Jpp master_rocky-44-g75b7c4f75
the software that should make you happy
Loading...
Searching...
No Matches
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

◆ main()

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>();
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);
93 zap['@'] = make_field(parameters) = JPARSER::initialised();
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;
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 if (inputFile.size() > 1u) { // sort input files
174
176
177 for (const string& file_name : inputFile) {
178
179 STATUS(file_name << '\r'); DEBUG(endl);
180
181 for (JSingleFileScanner<JEvent> in(file_name, 1); in.hasNext(); ) {
182
183 const JEvent* evt = in.next();
184
185 if (detector.getID() != evt->getDetectorID()) {
186 FATAL("Invalid detector identifier " << evt->getDetectorID() << " != " << evt->getDetectorID() << endl);
187 }
188
189 if (!evt->empty()) {
190 zmap[evt->begin()->getToE()] = file_name;
191 }
192 }
193 }
194 STATUS(endl);
195
196 inputFile.clear();
197
198 for (map<double, string>::const_iterator i = zmap.begin(); i != zmap.end(); ++i) {
199 inputFile.push_back(i->second);
200 }
201 }
202
203 outputFile.open();
204
205 outputFile.put(JMeta(argc, argv));
206 outputFile.put(parameters);
207
208 try { // limit RAM usage
209 dynamic_cast<JTreeWriterObjectOutput<JEvent>&>(*outputFile).getTreeWriter().SetAutoFlush(5000);
210 }
211 catch(const exception&) {}
212
213 int counter[] = { 0, 0 };
214
215 typedef deque<JEvent> buffer_type;
216
217 for (buffer_type zbuf; inputFile.hasNext(); ) {
218
219 STATUS(inputFile.getFilename() << '\r'); DEBUG(endl);
220
221 // read one file at a time
222
223 for (const string file_name = inputFile.getFilename(); inputFile.hasNext() && file_name == inputFile.getFilename(); ) {
224
225 const JEvent* evt = inputFile.next();
226
227 if (!evt->empty() && emitters.has(evt->getID())) {
228 if (utc(evt->begin()->getToA()) && utc(evt->rbegin()->getToA())) {
229 zbuf.push_back(*evt);
230 }
231 }
232 }
233
234 sort(zbuf.begin(), zbuf.end()); // sort according first time-of-emission
235
236 for (buffer_type::iterator p = zbuf.begin(), q; p != zbuf.end(); p = q) {
237
238 for (q = p; ++q != zbuf.end() && q->begin()->getToE() <= p->rbegin()->getToE() + parameters.Tmax_s; ) {}
239
240 if (q == zbuf.end()) {
241
242 if (inputFile.hasNext()) {
243
244 zbuf.erase(zbuf.begin(), p); // remove processed data and continue reading
245
246 break;
247 }
248 }
249
250 JEvent::overlap(p, q, parameters.deadTime_s); // empty overlapping events
251
252 if (getNumberOfEmitters(p,q) >= parameters.Nmin) {
253
254 h1.Fill((double) getNumberOfEmitters(p,q));
255
256 const JWeight getWeight(p, q);
257
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 * (parameters.Qmin <= 1.0 ? 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 }
284 }
285 }
286
287 if (getMinimumNumberOfEmitters(data.begin(), data.end()) >= parameters.Nmin) {
288
289 for (data_type::const_iterator hit = data.begin(); hit != data.end(); ++hit) {
290 HA[hit->getID()]->Fill(geometry.getIndex(hit->getString()), hit->getFloor(), 1.0);
291 }
292
293 // fit
294
295 const auto result = katoomba(data.begin(), data.end());
296
297 for (data_type::const_iterator hit = result.begin; hit != result.end; ++hit) {
298 HB[hit->getID()]->Fill(geometry.getIndex(hit->getString()), hit->getFloor(), 1.0);
299 }
300
301 if (debug >= debug_t) {
302
303 cout << "result:" << ' '
304 << FIXED(12,3) << result.chi2 << '/' << FIXED(6,1) << result.ndf << endl
305 << result.value << endl;
306
307 for (data_type::const_iterator hit = result.begin; hit != result.end; ++hit) {
308 cout << "hit: " << *hit << ' ' << FIXED(9,3) << katoomba.evaluator(result.value, *hit) << endl;
309 }
310 }
311
312 hn.Fill(log10(katoomba.gandalf.numberOfIterations));
313 h0.Fill(result.chi2 / result.ndf);
314
315 // store results
316
317 if ((parameters.chi2perNDF > 0.0 && result.chi2 / result.ndf <= +parameters.chi2perNDF) ||
318 (parameters.chi2perNDF < 0.0 && result.chi2 / result.ndf >= -parameters.chi2perNDF)) {
319
320 const JEvt evt = getEvt(JHead(detector.getID(),
321 result.getTimeRange(),
322 data .size(),
323 result.size(),
324 result.value.getN(),
325 result.ndf,
326 result.chi2),
327 result.value);
328
329 outputFile.put(evt);
330
331 if (selection.is_valid<JEvent>()) {
332
333 for (buffer_type::iterator i = p; i != q; ++i) {
334
335 JEvent out(*i); // event with fitted times of emission
336
337 const double toe = result.value.emission[JEKey(i->getID(), distance(p,i))].t1;
338
339 for (JEvent::iterator hit = out.begin(); hit != out.end(); ++hit) {
340 hit->setToE(toe);
341 }
342
343 outputFile.put(out);
344 }
345
346 }
347
348 counter[0] += 1;
349
350 } else {
351
352 WARNING(endl << "Event not written - chi2/NDF " << FIXED(12,3) << result.chi2 << '/' << FIXED(6,1) << result.ndf << endl);
353
354 counter[1] += 1;
355 }
356
357 } else {
358
359 WARNING(endl << "Event not accepted - minimum number of emitters " << getMinimumNumberOfEmitters(data.begin(), data.end()) << endl);
360
361 counter[1] += 1;
362 }
363 }
364 }
365 }
366 STATUS(endl);
367
368 STATUS("Number of events written / rejected: " << counter[0] << " / " << counter[1] << endl);
369
370 if (!squash) {
371
372 JMultipleFileScanner<JMeta> io(inputFile);
373
374 io >> outputFile;
375 }
376
377 outputFile.put(h0);
378 outputFile.put(h1);
379 outputFile.put(hn);
380
381 for (JManager<int, TH2D>::const_iterator i = HA.begin(); i != HA.end(); ++i) { outputFile.put(*(i->second)); }
382 for (JManager<int, TH2D>::const_iterator i = HB.begin(); i != HB.end(); ++i) { outputFile.put(*(i->second)); }
383
384 outputFile.close();
385}
string outputFile
#define DEBUG(A)
Message macros.
Definition JMessage.hh:62
#define STATUS(A)
Definition JMessage.hh:63
#define FATAL(A)
Definition JMessage.hh:67
int debug
debug level
Definition JSirene.cc:69
#define WARNING(A)
Definition JMessage.hh:65
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition JParser.hh:2142
#define MAKE_CSTRING(A)
Make C-string.
Definition JPrint.hh:72
std::vector< T >::difference_type distance(typename std::vector< T >::const_iterator first, typename PhysicsEvent::const_iterator< T > second)
Specialisation of STL distance.
Monte Carlo run header.
Definition JHead.hh:1236
Detector data structure.
Definition JDetector.hh:96
General exception.
Definition JException.hh:24
Utility class to parse command line options.
Definition JParser.hh:1698
Auxiliary class to manage set of compatible ROOT objects (e.g. histograms) using unique keys.
Definition JManager.hh:47
JTreeWriter< T > & getTreeWriter()
Get TreeWriter.
Object writing to file.
General purpose class for object reading from a list of file names.
Object reading from a list of files.
virtual bool hasNext() override
Check availability of next element.
bool has(const T &value) const
Test whether given value is present.
Range of values.
Definition JRange.hh:42
T getLength() const
Get length (difference between upper and lower limit).
Definition JRange.hh:289
T getLowerLimit() const
Get lower limit.
Definition JRange.hh:202
T getUpperLimit() const
Get upper limit.
Definition JRange.hh:213
size_t getMinimumNumberOfEmitters(T __begin, T __end)
Get minimum number of emitters for any string in data.
double getWeight(T __begin, T __end)
Get total weight of data points.
JContainer< std::vector< JTripod > > tripods_container
Definition JSydney.cc:78
JContainer< std::vector< JTransmitter > > transmitters_container
Definition JSydney.cc:80
JContainer< std::vector< JHydrophone > > hydrophones_container
Definition JSydney.cc:79
size_t getNumberOfEmitters(T __begin, T __end)
Get number of emitters.
static const JSoundVelocity getSoundVelocity(1541.0, -17.0e-3, -2000.0)
Function object for velocity of sound.
static JDetectorMechanics getMechanics
Function object to get string mechanics.
JEvt getEvt(const JHead &header, const JModel &model)
Get event.
floor_range getRangeOfFloors(const JDetector &detector)
Get range of floors.
void load(const std::string &file_name, JDetector &detector)
Load detector from input file.
static const JStringCounter getNumberOfStrings
Function object to count unique strings.
@ debug_t
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.
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
std::vector< JHitW0 > buffer_type
hits
Definition JPerth.cc:71
std::vector< event_type > data_type
Definition JPerth.cc:82
return result
Definition JPolint.hh:853
Auxiliary data structure for floating point format specification.
Definition JManip.hh:448
Detector file.
Definition JHead.hh:227
Emitter key.
Definition JEKey.hh:36
Acoustic emitter.
Definition JEmitter.hh:30
int getID() const
Get emitter identifier.
const int getDetectorID() const
Get detector identifier.
Acoustic event fit.
double Qmin
minimal quality transmission
double deadTime_s
dead time between events [s]
size_t Nmin
minimum number of emitters
double sigma_s
time-of-arrival resolution [s]
double Tmax_s
time window to combine events [s]
double chi2perNDF
maximal chi2/NDF to store event
Global fit of prameterised detector geometry to acoustics data.
Definition JGlobalfit.hh:29
Acoustics hit.
Template definition of fit function of acoustic model.
Implementation for depth dependend velocity of sound.
JSoundVelocity & set(const double z0)
Set depth.
Acoustic transmission identifier.
int getID() const
Get identifier.
Auxiliary data structure to unify weights of acoustics data according to the number of pings per emit...
Data structure for measured coincidence rates of all pairs of PMTs in optical module.
Definition JFitK40.hh:103
Auxiliary wrapper for I/O of container with optional comment (see JComment).
Definition JContainer.hh:42
Empty structure for specification of parser element that is initialised (i.e. does not require input)...
Definition JParser.hh:68
Auxiliary class for ROOT class selection.
bool is_valid() const
Get status of given data type.
Auxiliary class for defining the range of iterations of objects.
Definition JLimit.hh:45
static counter_type max()
Get maximum counter value.
Definition JLimit.hh:128
Auxiliary class for ROOT I/O of application specific meta data.
Definition JMeta.hh:72
General purpose class for hash map of unique keys.
Definition JHashMap.hh:75