Jpp 20.0.0-195-g190c9e876
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 JDetectorMechanics mechanics; // mechanical model data
79 JFitParameters parameters; // fit parameters
80 disable_container disable; // disable tansmissions
81 JROOTClassSelection selection = getROOTClassSelection<JTYPELIST<JEvent>::typelist>();
83 bool squash;
84 Long64_t autoflush;
85 int debug;
86
87 try {
88
89 JParser<> zap("Application to fit position calibration model to acoustic data.");
90
91 zap['f'] = make_field(inputFile, "output of JAcousticEventBuilder[.sh]");
92 zap['o'] = make_field(outputFile);
93 zap['n'] = make_field(numberOfEvents) = JLimit::max();
94 zap['a'] = make_field(detectorFile);
95 zap['@'] = make_field(parameters) = JPARSER::initialised();
96 zap['V'] = make_field(V, "sound velocity") = JPARSER::initialised();
97 zap['T'] = make_field(tripods, "tripod data");
98 zap['Y'] = make_field(transmitters, "transmitter data") = JPARSER::initialised();
99 zap['H'] = make_field(hydrophones, "hydrophone data") = JPARSER::initialised();
100 zap['M'] = make_field(mechanics, "mechanics data") = JPARSER::initialised();
101 zap['!'] = make_field(disable, "disable transmission") = JPARSER::initialised();
102 zap['C'] = make_field(selection,
103 "Precede name of data structure by a '+' or '-' "
104 "to add or remove data types in the output, respectively."
105 "\nROOT wildcards are accepted.") = JPARSER::initialised();
106 zap['u'] = make_field(utc, "select UTC time range") = JRange<double>();
107 zap['q'] = make_field(squash, "squash meta data");
108 zap['R'] = make_field(autoflush) = 5000;
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
119
120 try {
121 load(detectorFile, detector);
122 }
123 catch(const JException& error) {
124 FATAL(error);
125 }
126
128
129 JHashMap<int, JLocation> receivers;
131
132 for (JDetector::const_iterator i = detector.begin(); i != detector.end(); ++i) {
133 receivers[i->getID()] = i->getLocation();
134 }
135
136 for (tripods_container::const_iterator i = tripods.begin(); i != tripods.end(); ++i) {
137 emitters[i->getID()] = JEmitter(i->getID(), i->getUTMPosition() - detector.getUTMPosition());
138 }
139
140 for (transmitters_container::const_iterator i = transmitters.begin(); i != transmitters.end(); ++i) {
141 try {
142 emitters[i->getID()] = JEmitter(i->getID(), i->getPosition() + detector.getModule(i->getLocation()).getPosition());
143 }
144 catch(const exception&) {} // if no module available, discard transmitter
145 }
146
147 V.set(detector.getUTMZ()); // sound velocity at detector depth
148
149 JGeometry geometry(detector, mechanics, hydrophones);
150 JGlobalfit 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 if (inputFile.size() > 1u) { // sort input files
177
179
180 for (const string& file_name : inputFile) {
181
182 STATUS(file_name << '\r'); DEBUG(endl);
183
184 for (JSingleFileScanner<JEvent> in(file_name, 1); in.hasNext(); ) {
185
186 const JEvent* evt = in.next();
187
188 if (detector.getID() != evt->getDetectorID()) {
189 FATAL("Invalid detector identifier " << evt->getDetectorID() << " != " << evt->getDetectorID() << endl);
190 }
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 outputFile.open();
207
208 outputFile.put(JMeta(argc, argv));
209 outputFile.put(parameters);
210 outputFile.put(mechanics);
211
212 try { // limit RAM usage
213 dynamic_cast<JTreeWriterObjectOutput<JEvent>&>(*outputFile).getTreeWriter().SetAutoFlush(autoflush);
214 dynamic_cast<JTreeWriterObjectOutput<JEvent>&>(*outputFile).getTreeWriter().SetCacheSize();
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
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 * (parameters.Qmin <= 1.0 ? 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(detector.getID(),
326 result.getTimeRange(),
327 data .size(),
328 result.size(),
329 result.value.getN(),
330 result.ndf,
331 result.chi2,
332 katoomba.gandalf.numberOfIterations),
333 result.value);
334
335 outputFile.put(evt);
336
337 if (selection.is_valid<JEvent>()) {
338
339 for (buffer_type::iterator i = p; i != q; ++i) {
340
341 JEvent out(*i); // event with fitted times of emission
342
343 const double toe = result.value.emission[JEKey(i->getID(), distance(p,i))].t1;
344
345 for (JEvent::iterator hit = out.begin(); hit != out.end(); ++hit) {
346 hit->setToE(toe);
347 }
348
349 if (!out.empty()) {
350 outputFile.put(out);
351 }
352 }
353 }
354
355 counter[0] += 1;
356
357 } else {
358
359 WARNING(endl << "Event not written - chi2/NDF " << FIXED(12,3) << result.chi2 << '/' << FIXED(6,1) << result.ndf << endl);
360
361 counter[1] += 1;
362 }
363
364 } else {
365
366 WARNING(endl << "Event not accepted - minimum number of emitters " << getMinimumNumberOfEmitters(data.begin(), data.end()) << endl);
367
368 counter[1] += 1;
369 }
370 }
371 }
372 }
373 STATUS(endl);
374
375 STATUS("Number of events written / rejected: " << counter[0] << " / " << counter[1] << endl);
376
377 if (!squash) {
378
379 JMultipleFileScanner<JMeta> io(inputFile);
380
381 io >> outputFile;
382 }
383
384 outputFile.put(h0);
385 outputFile.put(h1);
386 outputFile.put(hn);
387
388 for (JManager<int, TH2D>::const_iterator i = HA.begin(); i != HA.end(); ++i) { outputFile.put(*(i->second)); }
389 for (JManager<int, TH2D>::const_iterator i = HB.begin(); i != HB.end(); ++i) { outputFile.put(*(i->second)); }
390
391 outputFile.close();
392}
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:74
#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.
JContainer< std::vector< JTripod > > tripods_container
Definition JSydney.cc:79
JContainer< std::vector< JTransmitter > > transmitters_container
Definition JSydney.cc:81
JContainer< std::vector< JHydrophone > > hydrophones_container
Definition JSydney.cc:80
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.
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:74
std::vector< event_type > data_type
Definition JPerth.cc:105
JRECONSTRUCTION::JWeight getWeight
return result
Definition JPolint.hh:862
Auxiliary data structure for floating point format specification.
Definition JManip.hh:448
Detector file.
Definition JHead.hh:227
Auxiliary data structure for mechanical model parameters with commented data.
Definition JMechanics.hh:38
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