Jpp master_rocky-44-g75b7c4f75
the software that should make you happy
Loading...
Searching...
No Matches
JPerth.cc File Reference

Program to determine inter-string time calibration. More...

#include <string>
#include <iostream>
#include <iomanip>
#include <type_traits>
#include <functional>
#include <future>
#include <mutex>
#include <thread>
#include <vector>
#include <set>
#include <memory>
#include <algorithm>
#include "km3net-dataformat/online/JDAQEvent.hh"
#include "km3net-dataformat/definitions/fitparameters.hh"
#include "JDAQ/JDAQEventIO.hh"
#include "JDetector/JDetector.hh"
#include "JDetector/JDetectorToolkit.hh"
#include "JDetector/JModuleRouter.hh"
#include "JDynamics/JDynamics.hh"
#include "JSupport/JMeta.hh"
#include "JSupport/JSupport.hh"
#include "JSupport/JMultipleFileScanner.hh"
#include "JSupport/JParallelFileScanner.hh"
#include "JSupport/JSummaryFileRouter.hh"
#include "JReconstruction/JHitW0.hh"
#include "JReconstruction/JEvt.hh"
#include "JReconstruction/JEvtToolkit.hh"
#include "JReconstruction/JMuonGandalfParameters_t.hh"
#include "JTrigger/JHitL0.hh"
#include "JTrigger/JBuildL0.hh"
#include "JFit/JFitToolkit.hh"
#include "JFit/JLine1Z.hh"
#include "JFit/JLine3Z.hh"
#include "JFit/JModel.hh"
#include "JFit/JGandalf.hh"
#include "JFit/JRegressorHelper.hh"
#include "JFit/JLine3ZRegressor.hh"
#include "JFit/JGradient.hh"
#include "JLang/JSTDObjectReader.hh"
#include "JLang/JVectorize.hh"
#include "JGeometry3D/JRotation3D.hh"
#include "JMath/JMath.hh"
#include "JMath/JQuantile_t.hh"
#include "JTools/JRange.hh"
#include "Jeep/JProperties.hh"
#include "Jeep/JParser.hh"
#include "Jeep/JMessage.hh"

Go to the source code of this file.

Classes

struct  JRECONSTRUCTION::event_type
 Auxiliary data structure to store data and fit in memory. More...
 
struct  JRECONSTRUCTION::JEditor
 Auxiliary class for editing time offset. More...
 
struct  JRECONSTRUCTION::JPerth
 Thread pool for fits to data. More...
 
struct  JRECONSTRUCTION::JPerth_t
 Auxiliary data structure for chi2 function object. More...
 

Namespaces

namespace  JRECONSTRUCTION
 Model fits to data.
 

Typedefs

typedef std::vector< JHitW0JRECONSTRUCTION::buffer_type
 hits
 
typedef std::map< int, buffer_typeJRECONSTRUCTION::map_type
 identifier -> hits
 
typedef std::vector< event_typeJRECONSTRUCTION::data_type
 
typedef JLANG::JSTDObjectReader< const event_typeJRECONSTRUCTION::input_type
 
typedef JFIT::JRegressorStorage< JFIT::JLine3Z, JFIT::JGandalfJRECONSTRUCTION::JRegressorStorage_t
 
typedef JFIT::JRegressor< JFIT::JLine3Z, JFIT::JGandalfJRECONSTRUCTION::JRegressor_t
 

Functions

int main (int argc, char **argv)
 

Detailed Description

Program to determine inter-string time calibration.

Author
mdejong

Definition in file JPerth.cc.

Function Documentation

◆ main()

int main ( int argc,
char ** argv )

Definition at line 281 of file JPerth.cc.

282{
283 using namespace std;
284 using namespace JPP;
285 using namespace KM3NETDAQ;
286
287 typedef JParallelFileScanner< JTypeList<JDAQEvent, JFIT::JEvt> > JParallelFileScanner_t;
288 typedef JParallelFileScanner_t::multi_pointer_type multi_pointer_type;
290 JACOUSTICS::JEvt>::typelist calibration_types;
291 typedef JMultipleFileScanner<calibration_types> JCalibration_t;
292
293 JMultipleFileScanner<> inputFile;
294 JLimit_t& numberOfEvents = inputFile.getLimit();
295 string detectorFile;
296 JCalibration_t calibrationFile;
297 double Tmax_s;
298 string pdfFile;
299 JMuonGandalfParameters_t parameters;
300 bool overwriteDetector;
301 set<int> strings; // string identifiers
302 set<int> modules; // module identifiers
303 int number_of_iterations = 1000;
304 int number_of_extra_steps = 0;
305 double epsilon = 1.0e-4;
306 double T_ns = 2.5; // step size
307 size_t threads; // number of threads
308 int debug;
309
310 const int DEFAULT_ID = -1;
311
312 try {
313
314 JProperties properties;
315
316 properties.insert(gmake_property(number_of_iterations));
317 properties.insert(gmake_property(number_of_extra_steps));
318 properties.insert(gmake_property(epsilon));
319 properties.insert(gmake_property(T_ns));
320
321 JParser<> zap("Program to determine inter-string time calibration.");
322
323 zap['f'] = make_field(inputFile);
324 zap['a'] = make_field(detectorFile);
325 zap['+'] = make_field(calibrationFile) = JPARSER::initialised();
326 zap['T'] = make_field(Tmax_s) = 100.0;
327 zap['n'] = make_field(numberOfEvents) = JLimit::max();
328 zap['P'] = make_field(pdfFile);
329 zap['@'] = make_field(parameters) = JPARSER::initialised();
330 zap['A'] = make_field(overwriteDetector, "overwrite detector file provided through '-a' with fitted time offsets.");
331 zap['S'] = make_field(strings, "string identifiers (" << DEFAULT_ID << " => all)") = JPARSER::initialised();
332 zap['M'] = make_field(modules, "module identifiers (" << DEFAULT_ID << " => all)") = JPARSER::initialised();
333 zap['%'] = make_field(properties, "fit parameters") = JPARSER::initialised();
334 zap['N'] = make_field(threads, "number of threads") = 1;
335 zap['d'] = make_field(debug) = 1;
336
337 zap(argc, argv);
338 }
339 catch(const exception& error) {
340 FATAL(error.what() << endl);
341 }
342
343
344 if (strings.empty() == modules.empty()) {
345 FATAL("Set either strings (option -S) or modules (option -M)." << endl);
346 }
347
349
350 try {
351 load(detectorFile, detector);
352 }
353 catch(const JException& error) {
354 FATAL(error);
355 }
356
357 unique_ptr<JDynamics> dynamics;
358
359 try {
360
361 dynamics.reset(new JDynamics(detector, Tmax_s));
362
363 dynamics->load(calibrationFile);
364 }
365 catch(const exception& error) {
366 if (!calibrationFile.empty()) {
367 FATAL(error.what());
368 }
369 }
370
371 const JModuleRouter router(dynamics ? dynamics->getDetector() : detector);
372
373 NOTICE("Reading PDFs... " << flush);
374
375 const JRegressorStorage_t storage(pdfFile, parameters.TTS_ns);
376
377 NOTICE("OK" << endl);
378
379 JRegressor_t::debug = debug;
380 JRegressor_t::T_ns.setRange(parameters.TMin_ns, parameters.TMax_ns);
381 JRegressor_t::Vmax_npe = parameters.VMax_npe;
382 JRegressor_t::MAXIMUM_ITERATIONS = parameters.NMax;
383
384
385 // preprocess input data
386
387 const JBuildL0<JHitL0> buildL0;
388
390
391 counter_type skip = inputFile.getLowerLimit();
392 counter_type counter = inputFile.getLowerLimit();
393
394 for (JMultipleFileScanner<>::const_iterator i = inputFile.begin(); i != inputFile.end(); ++i) {
395
396 JSummaryFileRouter summary(*i, parameters.R_Hz);
397
398 for (JParallelFileScanner_t in(*i); (skip -= in.skip(skip)) == 0 && in.hasNext() && counter != numberOfEvents; ++counter) {
399
400 STATUS("event: " << setw(10) << counter << '\r'); DEBUG(endl);
401
402 multi_pointer_type ps = in.next();
403
404 const JDAQEvent* tev = ps;
405 const JFIT::JEvt* evt = ps;
406
407 summary.update(*tev);
408
409 if (dynamics) {
410 dynamics->update(*tev);
411 }
412
413 if (!evt->empty()) {
414
415 vector<JHitL0> dataL0;
416
417 buildL0(*tev, router, true, back_inserter(dataL0));
418
419 const JFIT::JFit& fit = (*evt)[0];
420
421 const JRotation3D R (getDirection(fit));
422 const JLine1Z tz(getPosition (fit).rotate(R), fit.getT());
423 JRange<double> Z_m;
424
425 if (fit.hasW(JSTART_LENGTH_METRES) &&
426 fit.getW(JSTART_LENGTH_METRES) > 0.0) {
427 Z_m = JZRange(parameters.ZMin_m, parameters.ZMax_m + fit.getW(JSTART_LENGTH_METRES));
428 }
429
430 const JFIT::JModel<JLine1Z> match(tz, parameters.roadWidth_m, JRegressor_t::T_ns, Z_m);
431
432 // hit selection based on start value
433
434 buffer_type buffer;
435
436 for (vector<JHitL0>::const_iterator i = dataL0.begin(); i != dataL0.end(); ++i) {
437
438 JHitW0 hit(*i, summary.getRate(i->getPMTIdentifier()));
439
440 hit.rotate(R);
441
442 if (match(hit)) {
443 buffer.push_back(hit);
444 }
445 }
446
447 // select first hit
448
449 sort(buffer.begin(), buffer.end(), JHitL0::compare);
450
451 buffer_type::const_iterator __end = unique(buffer.begin(), buffer.end(), equal_to<JDAQPMTIdentifier>());
452
453 // re-organise data
454
455 map_type map;
456
457 for (buffer_type::const_iterator hit = buffer.begin(); hit != __end; ++hit) {
458
459 const JModule& module = router.getModule(hit->getModuleID());
460
461 if (!strings.empty()) { map[module.getString()].push_back(*hit); }
462 if (!modules.empty()) { map[module.getID()] .push_back(*hit); }
463 }
464
465 data.push_back({map, tz});
466 }
467 }
468 }
469 STATUS(endl);
470
471
472 // fit
473
474 JGradient fit(number_of_iterations, number_of_extra_steps, epsilon, 3);
475
476 if (strings.count(DEFAULT_ID) != 0) { strings = getStringIDs(detector); }
477 if (modules.count(DEFAULT_ID) != 0) { modules = getModuleIDs(detector); }
478
479 for (int id : strings) { fit.push_back(JModifier_t(MAKE_STRING("string: " << FILL(4,'0') << id << FILL()), new JEditor(data, id), T_ns)); }
480 for (int id : modules) { fit.push_back(JModifier_t(MAKE_STRING("module: " << setw(8) << id), new JEditor(data, id), T_ns)); }
481
482 const JPerth_t perth = {storage, data, threads};
483
484 const double chi2 = fit(perth);
485
486 STATUS("result: " << FIXED(12,6) << chi2 << ' ' << setw(6) << fit.numberOfIterations << endl);
487
488 for (size_t i = 0; i != fit.size(); ++i) {
489 {
490 const JEditor* p = dynamic_cast<const JEditor*>(fit[i].get());
491
492 if (p != NULL) { STATUS(fit[i].name << ' ' << FIXED(9,3) << p->t0 << " [ns]" << endl); }
493 }
494 }
495
496
497 if (overwriteDetector) {
498
499 detector.comment.add(JMeta(argc, argv));
500
502
503 for (size_t i = 0; i != fit.size(); ++i) {
504
505 const JEditor* p = dynamic_cast<const JEditor*>(fit[i].get());
506
507 if (p != NULL) {
508 calibration[p->id] = p->t0;
509 }
510 }
511
512 const double t0 = getAverage(get_values(calibration), 0.0);
513
514 for (JDetector::iterator module = detector.begin(); module != detector.end(); ++module) {
515
516 if (!module->empty()) {
517
518 map<int, double>::const_iterator p = (!strings.empty() ? calibration.find(module->getString()) :
519 !modules.empty() ? calibration.find(module->getID()) :
520 calibration.end());
521 if (p != calibration.end()) {
522 for (int pmt = 0; pmt != NUMBER_OF_PMTS; ++pmt) {
523 module->getPMT(pmt).addT0(p->second - t0);
524 }
525 }
526 }
527 }
528
529 NOTICE("Store calibration data on file " << detectorFile << endl);
530
531 store(detectorFile, detector);
532 }
533}
#define DEBUG(A)
Message macros.
Definition JMessage.hh:62
#define STATUS(A)
Definition JMessage.hh:63
#define NOTICE(A)
Definition JMessage.hh:64
#define FATAL(A)
Definition JMessage.hh:67
int debug
debug level
Definition JSirene.cc:69
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition JParser.hh:2142
#define MAKE_STRING(A)
Make string.
Definition JPrint.hh:63
#define gmake_property(A)
macros to convert (template) parameter to JPropertiesElement object
Detector data structure.
Definition JDetector.hh:96
int getString() const
Get string number.
Definition JLocation.hh:135
Router for direct addressing of module data in detector data structure.
Data structure for a composite optical module.
Definition JModule.hh:75
Utility class to parse parameter values.
Data structure for set of track fit results.
Data structure for track fit results with history and optional associated values.
const std::vector< double > & getW() const
Get associated values.
double getT() const
Get time.
bool hasW(const int i) const
Check availability of value.
Data structure for fit of straight line paralel to z-axis.
Definition JLine1Z.hh:29
General exception.
Definition JException.hh:24
int getID() const
Get identifier.
Definition JObjectID.hh:50
Utility class to parse command line options.
Definition JParser.hh:1698
Auxiliary class for a hit with background rate value.
Definition JHitW0.hh:23
General purpose class for object reading from a list of file names.
General purpose class for parallel reading of objects from a single file or multiple files.
File router for fast addressing of summary data.
Range of values.
Definition JRange.hh:42
Template L0 hit builder.
Definition JBuildL0.hh:38
const char * map
Definition elog.cc:87
static const int JSTART_LENGTH_METRES
distance between first and last hits in metres from JStart.cc
const double epsilon
JDirection3D getDirection(const Vec &dir)
Get direction.
JPosition3D getPosition(const Vec &pos)
Get position.
void load(const std::string &file_name, JDetector &detector)
Load detector from input file.
void store(const std::string &file_name, const JDetector &detector)
Store detector to output file.
std::set< int > getStringIDs(const JDetector &detector)
Get list of strings identifiers.
std::set< int > getModuleIDs(const JDetector &detector, const bool option=false)
Get list of modules identifiers.
JTOOLS::JRange< double > JZRange
const array_type< JValue_t > & get_values(const std::map< JKey_t, JValue_t, JComparator_t, JAllocator_t > &data)
Method to create array of values of map.
std::iterator_traits< T >::value_type getAverage(T __begin, T __end)
Get average.
Definition JMath.hh:494
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
Long64_t counter_type
Type definition for counter.
KM3NeT DAQ data structures and auxiliaries.
Definition DataQueue.cc:39
static const int NUMBER_OF_PMTS
Total number of PMTs in module.
Definition JDAQ.hh:26
Auxiliary data structure for sequence of same character.
Definition JManip.hh:330
Auxiliary data structure for floating point format specification.
Definition JManip.hh:448
Calibration.
Definition JHead.hh:330
Detector file.
Definition JHead.hh:227
Acoustic event fit.
Data structure for measured coincidence rates of all pairs of PMTs in optical module.
Definition JFitK40.hh:103
Orientation of module.
Dynamic detector calibration.
Definition JDynamics.hh:84
Conjugate gradient fit.
Definition JGradient.hh:75
Auxiliary class to match data points with given model.
Auxiliary data structure for editable parameter.
Definition JGradient.hh:49
Template data structure for storage for PDF tables.
Auxiliary class for recursive type list generation.
Definition JTypeList.hh:351
Empty structure for specification of parser element that is initialised (i.e. does not require input)...
Definition JParser.hh:68
Auxiliary class for editing time offset.
Definition JPerth.cc:91
double t0
time offset [ns]
Definition JPerth.cc:128
double TMin_ns
minimal time w.r.t. Cherenkov hypothesis [ns]
double TMax_ns
maximal time w.r.t. Cherenkov hypothesis [ns]
double VMax_npe
maximum number of of photo-electrons
Auxiliary data structure for chi2 function object.
Definition JPerth.cc:253
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
Auxiliary data structure for sorting of hits.
Definition JHitL0.hh:85