Jpp  18.0.0-rc.1
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 make a global fit of the detector geometry to acoustic data. More...

#include <string>
#include <iostream>
#include <iomanip>
#include <vector>
#include <set>
#include <map>
#include <deque>
#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 "JROOT/JROOTClassSelector.hh"
#include "JDetector/JDetector.hh"
#include "JDetector/JDetectorToolkit.hh"
#include "JDetector/JTripod.hh"
#include "JDetector/JTransmitter.hh"
#include "JDetector/JModule.hh"
#include "JDetector/JHydrophone.hh"
#include "JTools/JHashMap.hh"
#include "JTools/JRange.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_t.hh"
#include "JAcoustics/JEvent.hh"
#include "JAcoustics/JEvt.hh"
#include "JAcoustics/JEvtToolkit.hh"
#include "JAcoustics/JSupport.hh"
#include "JAcoustics/JTransmission_t.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 make a global fit of the detector geometry to acoustic data.


Author
mdejong

Definition in file software/JAcoustics/JKatoomba.cc.

Function Documentation

int main ( int  argc,
char **  argv 
)

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

60 {
61  using namespace std;
62  using namespace JPP;
63 
64  typedef JContainer< vector<JTripod> > tripods_container;
65  typedef JContainer< vector<JTransmitter> > transmitters_container;
66  typedef JContainer< vector<JHydrophone> > hydrophones_container;
67  typedef JContainer< set<JTransmission_t> > disable_container;
68 
71  string detectorFile;
72  JLimit_t& numberOfEvents = inputFile.getLimit();
73  JSoundVelocity V = getSoundVelocity; // default sound velocity
74  tripods_container tripods; // tripods
75  transmitters_container transmitters; // transmitters
76  hydrophones_container hydrophones; // hydrophones
77  JFitParameters parameters; // fit parameters
78  bool unify; // unify weighing of pings
79  disable_container disable; // disable tansmissions
80  JROOTClassSelection selection = getROOTClassSelection<JTYPELIST<JEvent>::typelist>();
81  double Tmax_s; // deadtime [s]
82  int debug;
83 
84  try {
85 
86  JParser<> zap("Application to fit position calibration model to acoustic data.");
87 
88  zap['f'] = make_field(inputFile, "output of JAcousticEventBuilder[.sh]");
89  zap['o'] = make_field(outputFile);
90  zap['n'] = make_field(numberOfEvents) = JLimit::max();
91  zap['a'] = make_field(detectorFile);
93  zap['V'] = make_field(V, "sound velocity") = JPARSER::initialised();
94  zap['T'] = make_field(tripods, "tripod data");
95  zap['Y'] = make_field(transmitters, "transmitter data") = JPARSER::initialised();
96  zap['H'] = make_field(hydrophones, "hydrophone data") = JPARSER::initialised();
97  zap['M'] = make_field(getMechanics, "mechanics data") = JPARSER::initialised();
98  zap['u'] = make_field(unify, "unify weighing of pings");
99  zap['!'] = make_field(disable, "disable transmission") = JPARSER::initialised();
100  zap['C'] = make_field(selection,
101  "Precede name of data structure by a '+' or '-' "
102  "to add or remove data types in the output, respectively."
103  "\nROOT wildcards are accepted.") = JPARSER::initialised();
104  zap['D'] = make_field(Tmax_s, "deadtime [s]") = 100.0e-3;
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 
115 
116  try {
117  load(detectorFile, detector);
118  }
119  catch(const JException& error) {
120  FATAL(error);
121  }
122 
124 
125  JHashMap<int, JLocation> receivers;
126  JHashMap<int, JEmitter> emitters;
127 
128  for (JDetector::const_iterator i = detector.begin(); i != detector.end(); ++i) {
129  receivers[i->getID()] = i->getLocation();
130  }
131 
132  for (tripods_container::const_iterator i = tripods.begin(); i != tripods.end(); ++i) {
133  emitters[i->getID()] = JEmitter(i->getID(),
134  i->getUTMPosition() - detector.getUTMPosition());
135  }
136 
137  for (transmitters_container::const_iterator i = transmitters.begin(); i != transmitters.end(); ++i) {
138  try {
139  emitters[i->getID()] = JEmitter(i->getID(),
140  i->getPosition() + detector.getModule(i->getLocation()).getPosition());
141  }
142  catch(const exception&) {
143  continue; // if no module available, discard transmitter
144  }
145  }
146 
147  V.set(detector.getUTMZ()); // sound velocity at detector depth
148 
149  JGeometry geometry(detector, hydrophones);
150  JKatoomba_t katoomba(geometry, V, parameters);
151 
154 
155  DEBUG(geometry);
156 
157  typedef vector<JHit> data_type;
158 
159  TH1D h0("chi2/NDF", NULL, 50000, 0.0, 1000.0);
160  TH1D h1("h1", NULL, 51, -0.5, 50.5);
161  TH1D hn("hn", NULL, 100, 0.0, 6.0);
162 
163  JManager<int, TH2D> HA(new TH2D("ha[%]", NULL,
165  range.getLength() + 1, range.getLowerLimit() - 0.5, range.getUpperLimit() + 0.5));
166 
167  JManager<int, TH2D> HB(new TH2D("hb[%]", NULL,
169  range.getLength() + 1, range.getLowerLimit() - 0.5, range.getUpperLimit() + 0.5));
170 
171  for (Int_t i = 1; i <= HA->GetXaxis()->GetNbins(); ++i) {
172  HA->GetXaxis()->SetBinLabel(i, MAKE_CSTRING(geometry.at(i-1).first));
173  HB->GetXaxis()->SetBinLabel(i, MAKE_CSTRING(geometry.at(i-1).first));
174  }
175 
176  string oid; // detector identifier
177 
178  { // sort input files
179  map<double, string> zmap;
180 
181  for (const string& file_name : inputFile) {
182 
183  STATUS(file_name << '\r'); DEBUG(endl);
184 
185  for (JMultipleFileScanner<JEvent> in(file_name, 1); in.hasNext(); ) {
186 
187  const JEvent* evt = in.next();
188 
189  if (oid == "")
190  oid = evt->getOID();
191  else if (oid != evt->getOID()) // consistency check
192  FATAL("Invalid detector identifier " << evt->getOID() << " != " << oid << endl);
193 
194  if (!evt->empty()) {
195  zmap[evt->begin()->getToE()] = file_name;
196  }
197  }
198  }
199  STATUS(endl);
200 
201  inputFile.clear();
202 
203  for (map<double, string>::const_iterator i = zmap.begin(); i != zmap.end(); ++i) {
204  inputFile.push_back(i->second);
205  }
206  }
207 
208  const JRange<double> unit(0.0, 1.0);
209 
210  outputFile.open();
211 
212  outputFile.put(JMeta(argc, argv));
213  outputFile.put(parameters);
214 
215  try { // limit RAM usage
216  dynamic_cast<JTreeWriterObjectOutput<JEvent>&>(*outputFile).getTreeWriter().SetAutoFlush(5000);
217  }
218  catch(const exception&) {}
219 
220  int counter[] = { 0, 0 };
221 
222  typedef deque<JEvent> buffer_type;
223 
224  for (buffer_type zbuf; inputFile.hasNext(); ) {
225 
226  STATUS(inputFile.getFilename() << '\r'); DEBUG(endl);
227 
228  // read one file at a time
229 
230  for (const string file_name = inputFile.getFilename(); inputFile.hasNext() && file_name == inputFile.getFilename(); ) {
231 
232  const JEvent* evt = inputFile.next();
233 
234  if (!evt->empty() && emitters.has(evt->getID())) {
235  zbuf.push_back(*evt);
236  }
237  }
238 
239  sort(zbuf.begin(), zbuf.end()); // sort according first time-of-emission
240 
241  for (buffer_type::iterator p = zbuf.begin(), q; p != zbuf.end(); p = q) {
242 
243  for (q = p; ++q != zbuf.end() && q->begin()->getToE() <= p->rbegin()->getToE() + parameters.Tmax_s; ) {}
244 
245  if (q == zbuf.end()) {
246 
247  if (inputFile.hasNext()) {
248 
249  zbuf.erase(zbuf.begin(), p); // remove processed data and continue reading
250 
251  break;
252  }
253  }
254 
255  overlap(p, q, Tmax_s); // empty overlapping events
256 
257  if (getNumberOfEmitters(p,q) >= parameters.Nmin) {
258 
259  h1.Fill((double) getNumberOfEmitters(p,q));
260 
261  map<int, int> numberOfPings;
262 
263  for (buffer_type::const_iterator i = p; i != q; ++i) {
264  numberOfPings[i->getID()] += 1;
265  }
266 
267  for (map<int, int>::const_iterator i = numberOfPings.begin(); i != numberOfPings.end(); ++i) {
268  DEBUG("Number of pings " << setw(2) << i->first << ' ' << setw(3) << i->second << endl);
269  }
270 
271  int minimum_number_of_pings = numeric_limits<int>::max();
272 
273  for (map<int, int>::const_iterator i = numberOfPings.begin(); i != numberOfPings.end(); ++i) {
274  minimum_number_of_pings = min(minimum_number_of_pings, i->second);
275  }
276 
277  set<int> buffer;
278  data_type data;
279 
280  for (buffer_type::iterator evt = p; evt != q; ++evt) {
281 
282  sort(evt->begin(), evt->end(), compare);
283 
284  JEvent::iterator __end = unique(evt->begin(), evt->end(), make_comparator(&JTransmission::getID, JComparison::eq()));
285 
286  const JEmitter& emitter = emitters[evt->getID()];
287  const double weight = (unify ? (double) minimum_number_of_pings / (double) numberOfPings[evt->getID()] : 1.0);
288 
289  for (JEvent::const_iterator i = evt->begin(); i != __end; ++i) {
290 
291  if (disable.count(JTransmission_t(evt->getID(), i->getID())) == 0 &&
292  disable.count(JTransmission_t(-1, i->getID())) == 0) {
293 
294  if (receivers.has(i->getID()) && geometry.hasLocation(receivers[i->getID()]) && i->getQ() >= parameters.Qmin * (unit(parameters.Qmin) ? i->getW() : 1.0)) {
295 
296  data.push_back(JHit(emitter,
297  distance(p,evt),
298  receivers[i->getID()],
299  i->getToA(),
300  parameters.sigma_s,
301  weight));
302 
303  buffer.insert(evt->getID());
304  }
305  }
306  }
307  }
308 
309  if (buffer.size() >= parameters.Nmin) {
310 
311  for (data_type::const_iterator hit = data.begin(); hit != data.end(); ++hit) {
312  HA[hit->getID()]->Fill(geometry.getIndex(hit->getString()), hit->getFloor(), 1.0);
313  }
314 
315  // fit
316 
317  const auto result = katoomba(data.begin(), data.end());
318 
319  for (data_type::const_iterator hit = result.begin; hit != result.end; ++hit) {
320  HB[hit->getID()]->Fill(geometry.getIndex(hit->getString()), hit->getFloor(), 1.0);
321  }
322 
323  if (debug >= debug_t) {
324 
325  cout << "result:" << ' '
326  << FIXED(12,3) << result.chi2 << '/' << FIXED(6,1) << result.ndf << endl
327  << result.value << endl;
328 
329  for (data_type::const_iterator hit = result.begin; hit != result.end; ++hit) {
330  cout << "hit: " << *hit << ' ' << FIXED(9,3) << katoomba.evaluator(result.value, *hit) << endl;
331  }
332  }
333 
334  hn.Fill(log10(katoomba.gandalf.numberOfIterations));
335  h0.Fill(result.chi2 / result.ndf);
336 
337  // store results
338 
339  if ((parameters.chi2perNDF > 0.0 && result.chi2 / result.ndf <= +parameters.chi2perNDF) ||
340  (parameters.chi2perNDF < 0.0 && result.chi2 / result.ndf >= -parameters.chi2perNDF)) {
341 
342  const JEvt evt = getEvt(JHead(oid,
343  result.getTimeRange(),
344  data .size(),
345  result.size(),
346  result.value.getN(),
347  result.ndf,
348  result.chi2),
349  result.value);
350 
351  outputFile.put(evt);
352 
353  if (selection.is_valid<JEvent>()) {
354 
355  for (buffer_type::iterator i = p; i != q; ++i) {
356 
357  JEvent out(*i); // event with fitted times of emission
358 
359  const double toe = result.value.emission[JEKey(i->getID(), distance(p,i))].t1;
360 
361  for (JEvent::iterator hit = out.begin(); hit != out.end(); ++hit) {
362  hit->setToE(toe);
363  }
364 
365  outputFile.put(out);
366  }
367 
368  }
369 
370  counter[0] += 1;
371 
372  } else {
373 
374  WARNING(endl << "Event not written: " << FIXED(12,3) << result.chi2 << '/' << FIXED(6,1) << result.ndf << endl);
375 
376  counter[1] += 1;
377  }
378  }
379  }
380  }
381  }
382  STATUS(endl);
383 
384  STATUS("Number of events written / rejected: " << counter[0] << " / " << counter[1] << endl);
385 
386  JMultipleFileScanner<JMeta> io(inputFile);
387 
388  io >> outputFile;
389 
390  outputFile.put(h0);
391  outputFile.put(h1);
392  outputFile.put(hn);
393 
394  for (JManager<int, TH2D>::const_iterator i = HA.begin(); i != HA.end(); ++i) { outputFile.put(*(i->second)); }
395  for (JManager<int, TH2D>::const_iterator i = HB.begin(); i != HB.end(); ++i) { outputFile.put(*(i->second)); }
396 
397  outputFile.close();
398 }
Worker class for complete fit procedure of acoustic model.
Definition: JKatoomba_t.hh:29
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:1514
General exception.
Definition: JException.hh:23
void overlap(T p, T q, const double Tmax_s)
Empty overlapping events.
debug
Definition: JMessage.hh:29
JComparator< JResult_t T::*, JComparison::lt > make_comparator(JResult_t T::*member)
Helper method to create comparator between values of data member.
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:72
#define STATUS(A)
Definition: JMessage.hh:63
Detector data structure.
Definition: JDetector.hh:89
Auxiliary class for ROOT class selection.
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:136
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
floor_range getRangeOfFloors(const JDetector &detector)
Get range of floors.
Acoustics hit.
size_t getNumberOfEmitters(T __begin, T __end)
Get number of emitters.
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.
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 event fit.
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:1989
set_variable E_E log10(E_{fit}/E_{#mu})"
const std::string & getOID() const
Get detector identifier.
return result
Definition: JPolint.hh:764
JPosition3D getPosition(const Vec &pos)
Get position.
static struct JACOUSTICS::@4 compare
Auxiliary data structure to sort transmissions.
Monte Carlo run header.
Definition: JHead.hh:1221
static const JStringCounter getNumberOfStrings
Function object to count unique strings.
Implementation for depth dependend velocity of sound.
#define FATAL(A)
Definition: JMessage.hh:67
z range($ZMAX-$ZMIN)< $MINIMAL_DZ." fi fi typeset -Z 4 STRING typeset -Z 2 FLOOR JPlot1D -f $
JTreeWriter< T > & getTreeWriter()
Get TreeWriter.
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.
Emitter key.
Definition: JEKey.hh:32
Acoustic event.
const JLimit & getLimit() const
Get limit.
Definition: JLimit.hh:84
do echo n Creating graphics for string $STRING for((FLOOR=$FIRST_FLOOR;$FLOOR<=$LAST_FLOOR;FLOOR+=1))
do set_variable DETECTOR_TXT $WORKDIR detector
int getID() const
Get identifier.
then fatal Wrong number of arguments fi set_variable DETECTOR $argv[1] set_variable INPUT_FILE $argv[2] eval JPrintDetector a $DETECTOR O IDENTIFIER eval JPrintDetector a $DETECTOR O SUMMARY JAcoustics sh $DETECTOR_ID source JAcousticsToolkit sh CHECK_EXIT_CODE typeset A EMITTERS get_tripods $WORKDIR tripod txt EMITTERS get_transmitters $WORKDIR transmitter txt EMITTERS for EMITTER in
Definition: JCanberra.sh:46
Acoustic transmission identifier.
then echo WARNING
Definition: JTuneHV.sh:91
Template definition of fit function of acoustic model.
Definition: JKatoomba.hh:105
int debug
debug level
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62