Jpp  17.3.2
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/JTransmitter.hh"
#include "JDetector/JHydrophone.hh"
#include "JDetector/JModule.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 60 of file examples/JAcoustics/JKatoomba.cc.

61 {
62  using namespace std;
63  using namespace JPP;
64 
65  typedef JContainer< vector<JTripod> > tripods_container;
66  typedef JContainer< vector<JTransmitter> > transmitters_container;
67  typedef JContainer< vector<JHydrophone> > hydrophones_container;
68  typedef JContainer< map<int, JWaveform> > waveform_container;
69  typedef JContainer< set<JTransmission_t> > disable_container;
70 
71  string detectorFile;
72  string outputFile;
73  int numberOfEvents;
74  map<int, int> numberOfPings;
75  JSoundVelocity V = getSoundVelocity; // default sound velocity
76  tripods_container tripods; // tripods
77  transmitters_container transmitters; // transmitters
78  hydrophones_container hydrophones; // hydrophones
79  JTimeRange T_s(-1.0e-4, +1.0e-4); // time range of background [s]
80  JFitParameters parameters; // fit parameters
81  int fit; // type of fit
82  bool unify; // unify weighing of pings
83  disable_container disable; // disable tansmissions
84  waveform_container waveforms; // emitter power and frequency
85  UInt_t seed;
86  int debug;
87 
88  numberOfPings[WILDCARD] = NUMBER_OF_PINGS;
89  waveforms [WILDCARD] = EMITTER_WAVEFORM;
90 
91  try {
92 
93  JParser<> zap("Example application to test fit of model to simulated acoustic data.");
94 
95  zap['a'] = make_field(detectorFile);
96  zap['o'] = make_field(outputFile);
97  zap['n'] = make_field(numberOfEvents);
99  zap['N'] = make_field(numberOfPings, "number of pings per emitter") = JPARSER::initialised();
100  zap['V'] = make_field(V, "sound velocity") = JPARSER::initialised();
101  zap['T'] = make_field(tripods, "tripod data");
102  zap['Y'] = make_field(transmitters, "transmitter data") = JPARSER::initialised();
103  zap['H'] = make_field(hydrophones, "hydrophone data") = JPARSER::initialised();
104  zap['F'] = make_field(fit, "fit type") = linear_t, simplex_t, gandalf_t;
105  zap['u'] = make_field(unify, "unify weighing of pings");
106  zap['!'] = make_field(disable, "disable transmission") = JPARSER::initialised();
107  zap['W'] = make_field(waveforms, "emitter power and frequency") = JPARSER::initialised();
108  zap['S'] = make_field(seed) = 0;
109  zap['d'] = make_field(debug) = 1;
110 
111  zap(argc, argv);
112  }
113  catch(const exception &error) {
114  FATAL(error.what() << endl);
115  }
116 
117 
118  gRandom->SetSeed(seed);
119 
120  if (numberOfEvents <= 0) { FATAL("Invalid number of events " << numberOfEvents << endl); }
121 
123 
124  try {
125  load(detectorFile, detector);
126  }
127  catch(const JException& error) {
128  FATAL(error);
129  }
130 
131  vector<JEmitter> emitters;
132 
133  for (tripods_container::const_iterator i = tripods.begin(); i != tripods.end(); ++i) {
134  emitters.push_back(JEmitter(i->getID(),
135  i->getUTMPosition() - detector.getUTMPosition()));
136  }
137 
138  for (transmitters_container::const_iterator i = transmitters.begin(); i != transmitters.end(); ++i) {
139  try {
140  emitters.push_back(JEmitter(i->getID(),
141  i->getPosition() + detector.getModule(i->getLocation()).getPosition()));
142  }
143  catch(const exception&) {
144  continue; // if no string available, discard transmitter
145  }
146  }
147 
148  if (detector.empty()) { FATAL("No modules in detector." << endl); }
149  if (emitters.empty()) { FATAL("No emitters in system." << endl); }
150 
151 
152  V.set(detector.getUTMZ()); // sound velocity at detector depth
153 
154  const double D_m = -detector.getUTMZ(); // depth [m]
155 
156 
157  const JGeometry geometry(detector, hydrophones);
158 
159  DEBUG(geometry);
160 
161  JKatoomba<JEstimator> estimator(geometry, V, parameters.option);
162  JKatoomba<JAbstractMinimiser> evaluator(geometry, V, parameters.option);
163  JKatoomba<JSimplex> simplex (geometry, V, parameters.option);
164  JKatoomba<JGandalf> gandalf (geometry, V, parameters.option);
165 
166  simplex.estimator.reset(getMEstimator(parameters.mestimator));
167  gandalf.estimator.reset(getMEstimator(parameters.mestimator));
168 
169  simplex.debug = debug;
170  gandalf.debug = debug;
171 
173 
174 
175  typedef JHit<JPDFGauss> hit_type;
176 
177 
178  TH1D h0("cpu", NULL, 100, 1.0, 7.0);
179  TH1D h1("chi2/NDF", NULL, 100, 0.0, 5.0);
180 
181  JManager<int, TH1D> H1(new TH1D("emitter[%]", NULL, 100, -5.0e-5, +5.0e-5));
182  JManager<int, TH2D> H2(new TH2D("string[%]", NULL, 100, -1.0, +1.0, 100, -1.0, +1.0));
183 
184 
185  JTimer timer;
186 
187  for (int number_of_events = 0, count = 0; number_of_events != numberOfEvents; ++number_of_events) {
188 
189  STATUS("event: " << setw(10) << number_of_events << '\r'); DEBUG(endl);
190 
191  JModel model; // randomised model data
192 
193  for (JDetector::const_iterator module = detector.begin(); module != detector.end(); ++module) {
194 
195  if (model.string[module->getString()] == JMODEL::JString()) {
196 
197  model.string[module->getString()] = JMODEL::JString(gRandom->Uniform(-2.0e-2, +2.0e-2),
198  gRandom->Uniform(-2.0e-2, +2.0e-2));
199  }
200  }
201 
202  DEBUG("Model" << endl << model << endl);
203 
204  // set module positions according actual model data
205 
206  for (JDetector::iterator module = detector.begin(); module != detector.end(); ++module) {
207  if (geometry.hasLocation(module->getLocation())) {
208  module->set(estimator.detector[module->getString()].getPosition(model.string[module->getString()], module->getFloor()));
209  }
210  }
211 
212 
213  // generate hits
214 
215  int minimum_number_of_pings = numeric_limits<int>::max();
216 
217  for (vector<JEmitter>::const_iterator emitter = emitters.begin(); emitter != emitters.end(); ++emitter) {
218  minimum_number_of_pings = min(minimum_number_of_pings, getValue(numberOfPings, emitter->getID(), numberOfPings[WILDCARD]));
219  }
220 
221  vector<hit_type> data;
222 
223  for (vector<JEmitter>::const_iterator emitter = emitters.begin(); emitter != emitters.end(); ++emitter) {
224 
225  const int number_of_pings = getValue(numberOfPings, emitter->getID(), numberOfPings[WILDCARD]);
226 
227  const double signal = (unify ? (double) minimum_number_of_pings / (double) number_of_pings : 1.0);
228 
229  for (int ping_counter = 0; ping_counter != number_of_pings; ++ping_counter) {
230 
231  ++count;
232 
233  const double toe_s = ping_counter * 10.0 + gRandom->Uniform(-1.0, +1.0);
234 
235  model.emitter[emitter->getID()][count] = toe_s;
236 
237  for (JDetector::const_iterator module = detector.begin(); module != detector.end(); ++module) {
238 
239  if (disable.count(JTransmission_t(emitter->getID(), module->getID())) == 0) {
240 
241  if (geometry.hasLocation(module->getLocation())) {
242 
243  const JWaveform& waveform = getValue(waveforms, emitter->getID(), waveforms[WILDCARD]);
244 
245  const double d_m = getDistance(module->getPosition(), emitter->getPosition());
246  const double toa_s = toe_s + V.getTime(d_m, emitter->getZ(), module->getZ());
247  const double Q = waveform.getQ(D_m, d_m);
248 
249  if (Q >= parameters.Qmin) {
250 
251  double t1_s = toa_s;
252 
253  if (gRandom->Rndm() >= parameters.background)
254  t1_s = gRandom->Gaus(toa_s, parameters.sigma_s);
255  else
256  t1_s = gRandom->Uniform(toa_s + T_s.getLowerLimit(),
257  toa_s + T_s.getUpperLimit());
258 
259  const hit_type hit(*emitter,
260  count,
261  module->getLocation(),
262  JPDFGauss(t1_s, parameters.sigma_s, signal, parameters.background));
263 
264  DEBUG("hit: " << hit << ' ' << FIXED(7,1) << Q << endl);
265 
266  data.push_back(hit);
267  }
268  }
269  }
270  }
271  }
272  }
273 
274 
275  timer.reset();
276  timer.start();
277 
278  JModel result = estimator(data.begin(), data.end()); // pre-fit
279  double chi2 = numeric_limits<double>::max();
280 
281  switch (fit) {
282 
283  case linear_t:
284 
285  chi2 = evaluator(result, data.begin(), data.end());
286  break;
287 
288  case simplex_t:
289 
290  simplex.value = result; // start value
291 
292  chi2 = simplex(data.begin(), data.end()) / simplex.estimator->getRho(1.0);
293  result = simplex.value;
294  break;
295 
296  case gandalf_t:
297 
298  gandalf.value = result; // start value
299 
300  chi2 = gandalf(data.begin(), data.end()) / gandalf.estimator->getRho(1.0);
301  result = gandalf.value;
302  break;
303 
304  default:
305  break;
306  }
307 
308  timer.stop();
309 
310  double W = 0.0;
311 
312  for (vector<hit_type>::const_iterator hit = data.begin(); hit != data.end(); ++hit) {
313  W += hit->getWeight();
314  }
315 
316  const int ndf = data.size() - result.getN();
317 
318  DEBUG("Final values" << endl
319  << FIXED(9,3) << chi2 << '/' << ndf << endl
320  << result << endl);
321 
322 
323  h0.Fill(log10((double) timer.usec_wall));
324  h1.Fill(chi2 / (double) (W - result.getN()));
325 
326  for (JHashMap<int, JMODEL::JString>::const_iterator i = model.string.begin(); i != model.string.end(); ++i) {
327  H2[i->first]->Fill((i->second.tx - result.string [i->first].tx) * 1.0e3, // mrad
328  (i->second.ty - result.string [i->first].ty) * 1.0e3); // mrad
329  }
330 
331  for (JHashMap<JEKey, JMODEL::JEmitter>::const_iterator i = model.emitter.begin(); i != model.emitter.end(); ++i) {
332  H1[i->first.getID()]->Fill(i->second.t1 - result.emitter[i->first].t1);
333  }
334  }
335  STATUS(endl);
336 
337 
338  TFile out(outputFile.c_str(), "recreate");
339 
340  out << h0
341  << h1
342  << H1
343  << H2;
344 
345  out.Write();
346  out.Close();
347 }
Utility class to parse command line options.
Definition: JParser.hh:1517
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:333
Template specialisation of fit function of acoustic model based on JSimplex minimiser.
Definition: JKatoomba.hh:462
*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 fatal Number of tripods
Definition: JFootprint.sh:45
Empty structure for specification of parser element that is initialised (i.e. does not require input)...
Definition: JParser.hh:83
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:106
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:553
Auxiliary class to manage set of compatible ROOT objects (e.g. histograms) using unique keys...
Definition: JManager.hh:43
Detector file.
Definition: JHead.hh:226
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:1993
set_variable E_E log10(E_{fit}/E_{#mu})"
return result
Definition: JPolint.hh:764
Auxiliary class for CPU timing and usage.
Definition: JTimer.hh:32
JPosition3D getPosition(const Vec &pos)
Get position.
double getQ(const double D_m, const double d_m) const
Get quality at given distance.
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.
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:272
then if[[!-f $DETECTOR]] then JDetector sh $DETECTOR fi cat $WORKDIR trigger_parameters txt<< EOFtrigger3DMuon.enabled=1;trigger3DMuon.numberOfHits=5;trigger3DMuon.gridAngle_deg=1;ctMin=0.0;TMaxLocal_ns=15.0;EOF set_variable TRIGGEREFFICIENCY_TRIGGERED_EVENTS_ONLY INPUT_FILES=() for((i=1;$i<=$NUMBER_OF_RUNS;++i));do JSirene.sh $DETECTOR $JPP_DATA/genhen.km3net_wpd_V2_0.evt.gz $WORKDIR/sirene_ ${i}.root JTriggerEfficiency.sh $DETECTOR $DETECTOR $WORKDIR/sirene_ ${i}.root $WORKDIR/trigger_efficiency_ ${i}.root $WORKDIR/trigger_parameters.txt $JPP_DATA/PMT_parameters.txt INPUT_FILES+=($WORKDIR/trigger_efficiency_ ${i}.root) done for ANGLE_DEG in $ANGLES_DEG[*];do set_variable SIGMA_NS 3.0 set_variable OUTLIERS 3 set_variable OUTPUT_FILE $WORKDIR/matrix\[${ANGLE_DEG}\deg\].root $JPP_DIR/examples/JReconstruction-f"$INPUT_FILES[*]"-o $OUTPUT_FILE-S ${SIGMA_NS}-A ${ANGLE_DEG}-O ${OUTLIERS}-d ${DEBUG}--!fiif[[$OPTION=="plot"]];then if((0));then for H1 in h0 h1;do JPlot1D-f"$WORKDIR/matrix["${^ANGLES_DEG}" deg].root:${H1}"-y"1 2e3"-Y-L TR-T""-\^"number of events [a.u.]"-> o chi2
Definition: JMatrixNZ.sh:106
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:203
do echo n Creating graphics for string $STRING for((FLOOR=$FIRST_FLOOR;$FLOOR<=$LAST_FLOOR;FLOOR+=1))
do set_variable DETECTOR_TXT $WORKDIR detector
Acoustic transmission identifier.
Utility class for emitter power and frequency.
Simplex fit.
Definition: JKatoomba.hh:107
int debug
debug level
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62
Gandalf fit.
Definition: JKatoomba.hh:108