Jpp  18.3.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  int debug;
82 
83  try {
84 
85  JParser<> zap("Application to fit position calibration model to acoustic data.");
86 
87  zap['f'] = make_field(inputFile, "output of JAcousticEventBuilder[.sh]");
88  zap['o'] = make_field(outputFile);
89  zap['n'] = make_field(numberOfEvents) = JLimit::max();
90  zap['a'] = make_field(detectorFile);
92  zap['V'] = make_field(V, "sound velocity") = JPARSER::initialised();
93  zap['T'] = make_field(tripods, "tripod data");
94  zap['Y'] = make_field(transmitters, "transmitter data") = JPARSER::initialised();
95  zap['H'] = make_field(hydrophones, "hydrophone data") = JPARSER::initialised();
96  zap['M'] = make_field(getMechanics, "mechanics data") = JPARSER::initialised();
97  zap['!'] = make_field(disable, "disable transmission") = JPARSER::initialised();
98  zap['C'] = make_field(selection,
99  "Precede name of data structure by a '+' or '-' "
100  "to add or remove data types in the output, respectively."
101  "\nROOT wildcards are accepted.") = JPARSER::initialised();
102  zap['d'] = make_field(debug) = 1;
103 
104  zap(argc, argv);
105  }
106  catch(const exception &error) {
107  FATAL(error.what() << endl);
108  }
109 
110 
112 
113  try {
114  load(detectorFile, detector);
115  }
116  catch(const JException& error) {
117  FATAL(error);
118  }
119 
121 
122  JHashMap<int, JLocation> receivers;
123  JHashMap<int, JEmitter> emitters;
124 
125  for (JDetector::const_iterator i = detector.begin(); i != detector.end(); ++i) {
126  receivers[i->getID()] = i->getLocation();
127  }
128 
129  for (tripods_container::const_iterator i = tripods.begin(); i != tripods.end(); ++i) {
130  emitters[i->getID()] = JEmitter(i->getID(), i->getUTMPosition() - detector.getUTMPosition());
131  }
132 
133  for (transmitters_container::const_iterator i = transmitters.begin(); i != transmitters.end(); ++i) {
134  try {
135  emitters[i->getID()] = JEmitter(i->getID(), i->getPosition() + detector.getModule(i->getLocation()).getPosition());
136  }
137  catch(const exception&) {} // if no module available, discard transmitter
138  }
139 
140  V.set(detector.getUTMZ()); // sound velocity at detector depth
141 
142  JGeometry geometry(detector, hydrophones);
143  JGlobalfit katoomba(geometry, V, parameters);
144 
147 
148  DEBUG(geometry);
149 
150  typedef vector<JHit> data_type;
151 
152  TH1D h0("chi2/NDF", NULL, 50000, 0.0, 1000.0);
153  TH1D h1("h1", NULL, 51, -0.5, 50.5);
154  TH1D hn("hn", NULL, 100, 0.0, 6.0);
155 
156  JManager<int, TH2D> HA(new TH2D("ha[%]", NULL,
158  range.getLength() + 1, range.getLowerLimit() - 0.5, range.getUpperLimit() + 0.5));
159 
160  JManager<int, TH2D> HB(new TH2D("hb[%]", NULL,
162  range.getLength() + 1, range.getLowerLimit() - 0.5, range.getUpperLimit() + 0.5));
163 
164  for (Int_t i = 1; i <= HA->GetXaxis()->GetNbins(); ++i) {
165  HA->GetXaxis()->SetBinLabel(i, MAKE_CSTRING(geometry.at(i-1).first));
166  HB->GetXaxis()->SetBinLabel(i, MAKE_CSTRING(geometry.at(i-1).first));
167  }
168 
169  string oid; // detector identifier
170 
171  if (inputFile.size() > 1u) { // sort input files
172 
173  map<double, string> zmap;
174 
175  for (const string& file_name : inputFile) {
176 
177  STATUS(file_name << '\r'); DEBUG(endl);
178 
179  for (JSingleFileScanner<JEvent> in(file_name, 1); in.hasNext(); ) {
180 
181  const JEvent* evt = in.next();
182 
183  if (oid == "")
184  oid = evt->getOID();
185  else if (oid != evt->getOID()) // consistency check
186  FATAL("Invalid detector identifier " << evt->getOID() << " != " << oid << endl);
187 
188  if (!evt->empty()) {
189  zmap[evt->begin()->getToE()] = file_name;
190  }
191  }
192  }
193  STATUS(endl);
194 
195  inputFile.clear();
196 
197  for (map<double, string>::const_iterator i = zmap.begin(); i != zmap.end(); ++i) {
198  inputFile.push_back(i->second);
199  }
200  }
201 
202  const JRange<double> unit(0.0, 1.0);
203 
204  outputFile.open();
205 
206  outputFile.put(JMeta(argc, argv));
207  outputFile.put(parameters);
208 
209  try { // limit RAM usage
210  dynamic_cast<JTreeWriterObjectOutput<JEvent>&>(*outputFile).getTreeWriter().SetAutoFlush(5000);
211  }
212  catch(const exception&) {}
213 
214  int counter[] = { 0, 0 };
215 
216  typedef deque<JEvent> buffer_type;
217 
218  for (buffer_type zbuf; inputFile.hasNext(); ) {
219 
220  STATUS(inputFile.getFilename() << '\r'); DEBUG(endl);
221 
222  // read one file at a time
223 
224  for (const string file_name = inputFile.getFilename(); inputFile.hasNext() && file_name == inputFile.getFilename(); ) {
225 
226  const JEvent* evt = inputFile.next();
227 
228  if (!evt->empty() && emitters.has(evt->getID())) {
229  zbuf.push_back(*evt);
230  }
231  }
232 
233  sort(zbuf.begin(), zbuf.end()); // sort according first time-of-emission
234 
235  for (buffer_type::iterator p = zbuf.begin(), q; p != zbuf.end(); p = q) {
236 
237  for (q = p; ++q != zbuf.end() && q->begin()->getToE() <= p->rbegin()->getToE() + parameters.Tmax_s; ) {}
238 
239  if (q == zbuf.end()) {
240 
241  if (inputFile.hasNext()) {
242 
243  zbuf.erase(zbuf.begin(), p); // remove processed data and continue reading
244 
245  break;
246  }
247  }
248 
249  JEvent::overlap(p, q, parameters.deadTime_s); // empty overlapping events
250 
251  if (getNumberOfEmitters(p,q) >= parameters.Nmin) {
252 
253  h1.Fill((double) getNumberOfEmitters(p,q));
254 
255  const JWeight getWeight(p, q);
256 
257  set<int> buffer;
258  data_type data;
259 
260  for (buffer_type::iterator evt = p; evt != q; ++evt) {
261 
262  sort(evt->begin(), evt->end(), JKatoomba<>::compare);
263 
264  JEvent::iterator __end = unique(evt->begin(), evt->end(), make_comparator(&JTransmission::getID, JComparison::eq()));
265 
266  const JEmitter& emitter = emitters [evt->getID()];
267  const double weight = getWeight(evt->getID());
268 
269  for (JEvent::const_iterator i = evt->begin(); i != __end; ++i) {
270 
271  if (disable.count(JTransmission_t(evt->getID(), i->getID())) == 0 &&
272  disable.count(JTransmission_t(-1, i->getID())) == 0) {
273 
274  if (receivers.has(i->getID()) && geometry.hasLocation(receivers[i->getID()]) && i->getQ() >= parameters.Qmin * (unit(parameters.Qmin) ? i->getW() : 1.0)) {
275 
276  data.push_back(JHit(emitter,
277  distance(p,evt),
278  receivers[i->getID()],
279  i->getToA(),
280  parameters.sigma_s,
281  weight));
282 
283  buffer.insert(evt->getID());
284  }
285  }
286  }
287  }
288 
289  if (buffer.size() >= parameters.Nmin) {
290 
291  for (data_type::const_iterator hit = data.begin(); hit != data.end(); ++hit) {
292  HA[hit->getID()]->Fill(geometry.getIndex(hit->getString()), hit->getFloor(), 1.0);
293  }
294 
295  // fit
296 
297  const auto result = katoomba(data.begin(), data.end());
298 
299  for (data_type::const_iterator hit = result.begin; hit != result.end; ++hit) {
300  HB[hit->getID()]->Fill(geometry.getIndex(hit->getString()), hit->getFloor(), 1.0);
301  }
302 
303  if (debug >= debug_t) {
304 
305  cout << "result:" << ' '
306  << FIXED(12,3) << result.chi2 << '/' << FIXED(6,1) << result.ndf << endl
307  << result.value << endl;
308 
309  for (data_type::const_iterator hit = result.begin; hit != result.end; ++hit) {
310  cout << "hit: " << *hit << ' ' << FIXED(9,3) << katoomba.evaluator(result.value, *hit) << endl;
311  }
312  }
313 
314  hn.Fill(log10(katoomba.gandalf.numberOfIterations));
315  h0.Fill(result.chi2 / result.ndf);
316 
317  // store results
318 
319  if ((parameters.chi2perNDF > 0.0 && result.chi2 / result.ndf <= +parameters.chi2perNDF) ||
320  (parameters.chi2perNDF < 0.0 && result.chi2 / result.ndf >= -parameters.chi2perNDF)) {
321 
322  const JEvt evt = getEvt(JHead(oid,
323  result.getTimeRange(),
324  data .size(),
325  result.size(),
326  result.value.getN(),
327  result.ndf,
328  result.chi2),
329  result.value);
330 
331  outputFile.put(evt);
332 
333  if (selection.is_valid<JEvent>()) {
334 
335  for (buffer_type::iterator i = p; i != q; ++i) {
336 
337  JEvent out(*i); // event with fitted times of emission
338 
339  const double toe = result.value.emission[JEKey(i->getID(), distance(p,i))].t1;
340 
341  for (JEvent::iterator hit = out.begin(); hit != out.end(); ++hit) {
342  hit->setToE(toe);
343  }
344 
345  outputFile.put(out);
346  }
347 
348  }
349 
350  counter[0] += 1;
351 
352  } else {
353 
354  WARNING(endl << "Event not written: " << FIXED(12,3) << result.chi2 << '/' << FIXED(6,1) << result.ndf << endl);
355 
356  counter[1] += 1;
357  }
358  }
359  }
360  }
361  }
362  STATUS(endl);
363 
364  STATUS("Number of events written / rejected: " << counter[0] << " / " << counter[1] << endl);
365 
366  JMultipleFileScanner<JMeta> io(inputFile);
367 
368  io >> outputFile;
369 
370  outputFile.put(h0);
371  outputFile.put(h1);
372  outputFile.put(hn);
373 
374  for (JManager<int, TH2D>::const_iterator i = HA.begin(); i != HA.end(); ++i) { outputFile.put(*(i->second)); }
375  for (JManager<int, TH2D>::const_iterator i = HB.begin(); i != HB.end(); ++i) { outputFile.put(*(i->second)); }
376 
377  outputFile.close();
378 }
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: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:78
double getWeight(T __begin, T __end)
Get total weight of data points.
Definition: JKatoomba_t.hh:61
JContainer< std::vector< JTransmitter > > transmitters_container
Definition: JSydney.cc:79
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: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.
std::vector< JHitW0 > buffer_type
hits
Definition: JPerth.cc:67
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:78
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:1989
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:77
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.
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
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:48
Acoustic transmission identifier.
Template definition of fit function of acoustic model.
Definition: JKatoomba_t.hh:146
int debug
debug level
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62