Jpp  18.6.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/JSingleFileScanner.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/JGlobalfit.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 60 of file software/JAcoustics/JKatoomba.cc.

61 {
62  using namespace std;
63  using namespace JPP;
64 
68  typedef JContainer< set<JTransmission_t> > disable_container;
69 
72  string detectorFile;
73  JLimit_t& numberOfEvents = inputFile.getLimit();
74  JSoundVelocity V = getSoundVelocity; // default sound velocity
75  tripods_container tripods; // tripods
76  transmitters_container transmitters; // transmitters
77  hydrophones_container hydrophones; // hydrophones
78  JFitParameters parameters; // fit parameters
79  disable_container disable; // disable tansmissions
80  JROOTClassSelection selection = getROOTClassSelection<JTYPELIST<JEvent>::typelist>();
81  JRange<double> utc;
82  bool squash;
83  int debug;
84 
85  try {
86 
87  JParser<> zap("Application to fit position calibration model to acoustic data.");
88 
89  zap['f'] = make_field(inputFile, "output of JAcousticEventBuilder[.sh]");
90  zap['o'] = make_field(outputFile);
91  zap['n'] = make_field(numberOfEvents) = JLimit::max();
92  zap['a'] = make_field(detectorFile);
94  zap['V'] = make_field(V, "sound velocity") = JPARSER::initialised();
95  zap['T'] = make_field(tripods, "tripod data");
96  zap['Y'] = make_field(transmitters, "transmitter data") = JPARSER::initialised();
97  zap['H'] = make_field(hydrophones, "hydrophone data") = JPARSER::initialised();
98  zap['M'] = make_field(getMechanics, "mechanics data") = JPARSER::initialised();
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['u'] = make_field(utc, "select UTC time range") = JRange<double>();
105  zap['q'] = make_field(squash, "squash meta data");
106  zap['d'] = make_field(debug) = 1;
107 
108  zap(argc, argv);
109  }
110  catch(const exception &error) {
111  FATAL(error.what() << endl);
112  }
113 
114 
116 
117  try {
118  load(detectorFile, detector);
119  }
120  catch(const JException& error) {
121  FATAL(error);
122  }
123 
125 
126  JHashMap<int, JLocation> receivers;
127  JHashMap<int, JEmitter> emitters;
128 
129  for (JDetector::const_iterator i = detector.begin(); i != detector.end(); ++i) {
130  receivers[i->getID()] = i->getLocation();
131  }
132 
133  for (tripods_container::const_iterator i = tripods.begin(); i != tripods.end(); ++i) {
134  emitters[i->getID()] = JEmitter(i->getID(), 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(), i->getPosition() + detector.getModule(i->getLocation()).getPosition());
140  }
141  catch(const exception&) {} // if no module available, discard transmitter
142  }
143 
144  V.set(detector.getUTMZ()); // sound velocity at detector depth
145 
146  JGeometry geometry(detector, hydrophones);
147  JGlobalfit katoomba(geometry, V, parameters);
148 
151 
152  DEBUG(geometry);
153 
154  typedef vector<JHit> data_type;
155 
156  TH1D h0("chi2/NDF", NULL, 50000, 0.0, 1000.0);
157  TH1D h1("h1", NULL, 51, -0.5, 50.5);
158  TH1D hn("hn", NULL, 100, 0.0, 6.0);
159 
160  JManager<int, TH2D> HA(new TH2D("ha[%]", NULL,
162  range.getLength() + 1, range.getLowerLimit() - 0.5, range.getUpperLimit() + 0.5));
163 
164  JManager<int, TH2D> HB(new TH2D("hb[%]", NULL,
166  range.getLength() + 1, range.getLowerLimit() - 0.5, range.getUpperLimit() + 0.5));
167 
168  for (Int_t i = 1; i <= HA->GetXaxis()->GetNbins(); ++i) {
169  HA->GetXaxis()->SetBinLabel(i, MAKE_CSTRING(geometry.at(i-1).first));
170  HB->GetXaxis()->SetBinLabel(i, MAKE_CSTRING(geometry.at(i-1).first));
171  }
172 
173  string oid; // detector identifier
174 
175  if (inputFile.size() > 1u) { // sort input files
176 
177  map<double, string> zmap;
178 
179  for (const string& file_name : inputFile) {
180 
181  STATUS(file_name << '\r'); DEBUG(endl);
182 
183  for (JSingleFileScanner<JEvent> in(file_name, 1); in.hasNext(); ) {
184 
185  const JEvent* evt = in.next();
186 
187  if (oid == "")
188  oid = evt->getOID();
189  else if (oid != evt->getOID()) // consistency check
190  FATAL("Invalid detector identifier " << evt->getOID() << " != " << oid << endl);
191 
192  if (!evt->empty()) {
193  zmap[evt->begin()->getToE()] = file_name;
194  }
195  }
196  }
197  STATUS(endl);
198 
199  inputFile.clear();
200 
201  for (map<double, string>::const_iterator i = zmap.begin(); i != zmap.end(); ++i) {
202  inputFile.push_back(i->second);
203  }
204  }
205 
206  const JRange<double> unit(0.0, 1.0);
207 
208  outputFile.open();
209 
210  outputFile.put(JMeta(argc, argv));
211  outputFile.put(parameters);
212 
213  try { // limit RAM usage
214  dynamic_cast<JTreeWriterObjectOutput<JEvent>&>(*outputFile).getTreeWriter().SetAutoFlush(5000);
215  }
216  catch(const exception&) {}
217 
218  int counter[] = { 0, 0 };
219 
220  typedef deque<JEvent> buffer_type;
221 
222  for (buffer_type zbuf; inputFile.hasNext(); ) {
223 
224  STATUS(inputFile.getFilename() << '\r'); DEBUG(endl);
225 
226  // read one file at a time
227 
228  for (const string file_name = inputFile.getFilename(); inputFile.hasNext() && file_name == inputFile.getFilename(); ) {
229 
230  const JEvent* evt = inputFile.next();
231 
232  if (!evt->empty() && emitters.has(evt->getID())) {
233  if (utc(evt->begin()->getToA()) && utc(evt->rbegin()->getToA())) {
234  zbuf.push_back(*evt);
235  }
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  JEvent::overlap(p, q, parameters.deadTime_s); // empty overlapping events
256 
257  if (getNumberOfEmitters(p,q) >= parameters.Nmin) {
258 
259  h1.Fill((double) getNumberOfEmitters(p,q));
260 
261  const JWeight getWeight(p, q);
262 
263  data_type data;
264 
265  for (buffer_type::iterator evt = p; evt != q; ++evt) {
266 
267  sort(evt->begin(), evt->end(), JKatoomba<>::compare);
268 
269  JEvent::iterator __end = unique(evt->begin(), evt->end(), make_comparator(&JTransmission::getID, JComparison::eq()));
270 
271  const JEmitter& emitter = emitters [evt->getID()];
272  const double weight = getWeight(evt->getID());
273 
274  for (JEvent::const_iterator i = evt->begin(); i != __end; ++i) {
275 
276  if (disable.count(JTransmission_t(evt->getID(), i->getID())) == 0 &&
277  disable.count(JTransmission_t(-1, i->getID())) == 0) {
278 
279  if (receivers.has(i->getID()) && geometry.hasLocation(receivers[i->getID()]) && i->getQ() >= parameters.Qmin * (unit(parameters.Qmin) ? i->getW() : 1.0)) {
280 
281  data.push_back(JHit(emitter,
282  distance(p,evt),
283  receivers[i->getID()],
284  i->getToA(),
285  parameters.sigma_s,
286  weight));
287  }
288  }
289  }
290  }
291 
292  if (getMinimumNumberOfEmitters(data.begin(), data.end()) >= parameters.Nmin) {
293 
294  for (data_type::const_iterator hit = data.begin(); hit != data.end(); ++hit) {
295  HA[hit->getID()]->Fill(geometry.getIndex(hit->getString()), hit->getFloor(), 1.0);
296  }
297 
298  // fit
299 
300  const auto result = katoomba(data.begin(), data.end());
301 
302  for (data_type::const_iterator hit = result.begin; hit != result.end; ++hit) {
303  HB[hit->getID()]->Fill(geometry.getIndex(hit->getString()), hit->getFloor(), 1.0);
304  }
305 
306  if (debug >= debug_t) {
307 
308  cout << "result:" << ' '
309  << FIXED(12,3) << result.chi2 << '/' << FIXED(6,1) << result.ndf << endl
310  << result.value << endl;
311 
312  for (data_type::const_iterator hit = result.begin; hit != result.end; ++hit) {
313  cout << "hit: " << *hit << ' ' << FIXED(9,3) << katoomba.evaluator(result.value, *hit) << endl;
314  }
315  }
316 
317  hn.Fill(log10(katoomba.gandalf.numberOfIterations));
318  h0.Fill(result.chi2 / result.ndf);
319 
320  // store results
321 
322  if ((parameters.chi2perNDF > 0.0 && result.chi2 / result.ndf <= +parameters.chi2perNDF) ||
323  (parameters.chi2perNDF < 0.0 && result.chi2 / result.ndf >= -parameters.chi2perNDF)) {
324 
325  const JEvt evt = getEvt(JHead(oid,
326  result.getTimeRange(),
327  data .size(),
328  result.size(),
329  result.value.getN(),
330  result.ndf,
331  result.chi2),
332  result.value);
333 
334  outputFile.put(evt);
335 
336  if (selection.is_valid<JEvent>()) {
337 
338  for (buffer_type::iterator i = p; i != q; ++i) {
339 
340  JEvent out(*i); // event with fitted times of emission
341 
342  const double toe = result.value.emission[JEKey(i->getID(), distance(p,i))].t1;
343 
344  for (JEvent::iterator hit = out.begin(); hit != out.end(); ++hit) {
345  hit->setToE(toe);
346  }
347 
348  outputFile.put(out);
349  }
350 
351  }
352 
353  counter[0] += 1;
354 
355  } else {
356 
357  WARNING(endl << "Event not written - chi2/NDF " << FIXED(12,3) << result.chi2 << '/' << FIXED(6,1) << result.ndf << endl);
358 
359  counter[1] += 1;
360  }
361 
362  } else {
363 
364  WARNING(endl << "Event not accepted - minimum number of emitters " << getMinimumNumberOfEmitters(data.begin(), data.end()) << endl);
365 
366  counter[1] += 1;
367  }
368  }
369  }
370  }
371  STATUS(endl);
372 
373  STATUS("Number of events written / rejected: " << counter[0] << " / " << counter[1] << endl);
374 
375  if (!squash) {
376 
377  JMultipleFileScanner<JMeta> io(inputFile);
378 
379  io >> outputFile;
380  }
381 
382  outputFile.put(h0);
383  outputFile.put(h1);
384  outputFile.put(hn);
385 
386  for (JManager<int, TH2D>::const_iterator i = HA.begin(); i != HA.end(); ++i) { outputFile.put(*(i->second)); }
387  for (JManager<int, TH2D>::const_iterator i = HB.begin(); i != HB.end(); ++i) { outputFile.put(*(i->second)); }
388 
389  outputFile.close();
390 }
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:1711
General exception.
Definition: JException.hh:24
#define WARNING(A)
Definition: JMessage.hh:65
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.
std::vector< event_type > data_type
Definition: JPerth.cc:81
double getWeight(T __begin, T __end)
Get total weight of data points.
Definition: JKatoomba_t.hh:61
size_t getMinimumNumberOfEmitters(T __begin, T __end)
Get minimum number of emitters for any string in data.
JContainer< std::vector< JTransmitter > > transmitters_container
Definition: JSydney.cc:80
static JDetectorMechanics getMechanics
Function object to get string mechanics.
Definition: JMechanics.hh:242
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
Auxiliary data structure to unify weights of acoustics data according to the number of pings per emit...
#define STATUS(A)
Definition: JMessage.hh:63
Detector data structure.
Definition: JDetector.hh:89
Auxiliary class for ROOT class selection.
*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:84
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.
std::vector< JHitW0 > buffer_type
hits
Definition: JPerth.cc:70
Acoustics hit.
Global fit of prameterised detector geometry to acoustics data.
Definition: JGlobalfit.hh:29
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
JContainer< std::vector< JHydrophone > > hydrophones_container
Definition: JSydney.cc:79
Acoustic event fit.
Monte Carlo run header.
Definition: JHead.hh:1234
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:2158
set_variable E_E log10(E_{fit}/E_{#mu})"
const std::string & getOID() const
Get detector identifier.
JPosition3D getPosition(const Vec &pos)
Get position.
JContainer< std::vector< JTripod > > tripods_container
Definition: JSydney.cc:78
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.
then fatal The output file must have the wildcard in the e g root fi 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:48
Data structure for measured coincidence rates of all pairs of PMTs in optical module.
Definition: JFitK40.hh:99
JEvt getEvt(const JHead &header, const JModel &model)
Get event.
Emitter key.
Definition: JEKey.hh:32
Acoustic event.
Object reading from a list of files.
const JLimit & getLimit() const
Get limit.
Definition: JLimit.hh:84
bool overlap(const JRange< T, JComparator_t > &first, const JRange< T, JComparator_t > &second)
Test overlap between ranges.
Definition: JRange.hh:641
do set_variable DETECTOR_TXT $WORKDIR detector
int getID() const
Get emitter identifier.
double u[N+1]
Definition: JPolint.hh:865
Acoustic transmission identifier.
Template definition of fit function of acoustic model.
Definition: JKatoomba_t.hh:77
int debug
debug level
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62