#include #include #include #include "km3net-dataformat/offline/Head.hh" #include "km3net-dataformat/offline/MultiHead.hh" #include "km3net-dataformat/online/JDAQ.hh" #include "km3net-dataformat/online/JDAQClock.hh" #include "JDAQ/JDAQTimesliceIO.hh" #include "JDAQ/JDAQEventIO.hh" #include "JDAQ/JDAQSummarysliceIO.hh" #include "JDetector/JDetector.hh" #include "JDetector/JDetectorToolkit.hh" #include "JDetector/JModuleRouter.hh" #include "JTrigger/JHit.hh" #include "JTrigger/JHitToolkit.hh" #include "JTrigger/JSuperFrame1D.hh" #include "JTrigger/JSuperFrame2D.hh" #include "JTrigger/JTimeslice.hh" #include "JTrigger/JBuildL1.hh" #include "JTrigger/JBuildL2.hh" #include "JTrigger/JTimesliceL1.hh" #include "JTrigger/JTriggerParameters.hh" #include "JTrigger/JTriggerToolkit.hh" #include "JLang/JObjectMultiplexer.hh" #include "JROOT/JROOTClassSelector.hh" #include "JSupport/JMonteCarloFileSupportkit.hh" #include "JSupport/JSingleFileScanner.hh" #include "JSupport/JFileRecorder.hh" #include "JSupport/JSupport.hh" #include "JSupport/JMeta.hh" #include "JSupport/JTriggerParametersSupportkit.hh" #include "Jeep/JParser.hh" #include "Jeep/JMessage.hh" /** * \file * Auxiliary program to re-process KM3NETDAQ::JDAQTimeslice data. * \author mdejong */ int main(int argc, char **argv) { using namespace std; using namespace JPP; using namespace KM3NETDAQ; typedef JAllTypes_t typelist; JSingleFileScanner inputFile; JFileRecorder outputFile; JLimit_t& numberOfEvents = inputFile.getLimit(); string detectorFile; JTriggerParameters parameters; bool reuse_parameters; JROOTClassSelector selector; int debug; try { JParser<> zap("Auxiliary program to re-process time slice data."); zap['f'] = make_field(inputFile); zap['o'] = make_field(outputFile) = "timeslice.root"; zap['n'] = make_field(numberOfEvents) = JLimit::max(); zap['a'] = make_field(detectorFile); zap['@'] = make_field(parameters) = JPARSER::initialised(); zap['U'] = make_field(reuse_parameters); zap['C'] = make_field(selector, "timeslice selector, e.g. JDAQTimesliceL1.") = getROOTClassSelection(); zap['d'] = make_field(debug) = 0; zap(argc, argv); } catch(const exception &error) { FATAL(error.what() << endl); } setDAQLongprint(debug >= JEEP::debug_t); DEBUG("Frame time [ms] " << getFrameTime() * 1e-6 << endl); DEBUG("Reset time [ms] " << getRTS() * 1e-6 << endl); DEBUG("Trigger" << endl << parameters << endl); JDetector detector; try { load(detectorFile, detector); } catch(const JException& error) { FATAL(error); } if (reuse_parameters) { try { parameters = getTriggerParameters(inputFile); NOTICE("Set trigger parameters from input." << endl); } catch(const JException& error) { FATAL("No trigger parameters from input." << endl); } } // detector if (parameters.disableHighRateVeto) { NOTICE("Disabling high-rate veto of all PMTs." << endl); detector.setPMTStatus(HIGH_RATE_VETO_DISABLE); } parameters.set(getMaximalDistance(detector)); const JModuleRouter moduleRouter(detector); if (parameters.writeSummary()) { WARNING("Discard writeSummary option during reprocesing of data." << endl); } //typedef JHit hit_type; //typedef int hit_type; typedef double hit_type; typedef JSuperFrame1D JSuperFrame1D_t; typedef JSuperFrame2D JSuperFrame2D_t; typedef JTimeslice JTimeslice_t; typedef JBuildL1 JBuildL1_t; typedef JBuildL2 JBuildL2_t; const JBuildL1_t buildL1(parameters); const JBuildL2_t buildL2(parameters.L2); const JBuildL2_t buildSN(parameters.SN); JTimesliceRouter timesliceRouter(parameters.numberOfBins); outputFile.open(); if (!outputFile.is_open()) { FATAL("Error opening file " << outputFile << endl); } outputFile.put(JMeta(argc, argv)); outputFile.put(parameters); JObjectMultiplexer in(inputFile, selector); counter_type counter = 0; for ( ; in.hasNext() && counter != inputFile.getLimit(); ++counter) { STATUS("event: " << setw(10) << counter << '\r'); DEBUG(endl); const JDAQTimeslice* timeslice = in.next(); DEBUG(*timeslice << endl); timesliceRouter.configure(*timeslice); JTimeslice_t timesliceL0(timeslice->getDAQChronometer()); JTimeslice_t timesliceL1(timeslice->getDAQChronometer()); JTimeslice_t timesliceL2(timeslice->getDAQChronometer()); JTimeslice_t timesliceSN(timeslice->getDAQChronometer()); for (JDAQTimeslice::const_iterator super_frame = timeslice->begin(); super_frame != timeslice->end(); ++super_frame) { if (moduleRouter.hasModule(super_frame->getModuleID())) { // calibration const JModule& module = moduleRouter.getModule(super_frame->getModuleID()); JSuperFrame2D_t& buffer = JSuperFrame2D_t::demultiplex(*super_frame, module); // Apply high-rate veto buffer.applyHighRateVeto(parameters.highRateVeto_Hz); // L0 timesliceL0.push_back(JSuperFrame1D_t(buffer)); // L1 timesliceL1.push_back(JSuperFrame1D_t(super_frame->getDAQChronometer(), super_frame->getModuleIdentifier(), module.getPosition())); buildL1(*timesliceL0.rbegin(), back_inserter(*timesliceL1.rbegin())); // L2 timesliceL2.push_back(JSuperFrame1D_t(super_frame->getDAQChronometer(), super_frame->getModuleIdentifier(), module.getPosition())); buildL2(buffer, *timesliceL1.rbegin(), back_inserter(*timesliceL2.rbegin())); // SN timesliceSN.push_back(JSuperFrame1D_t(super_frame->getDAQChronometer(), super_frame->getModuleIdentifier(), module.getPosition())); buildSN(buffer, *timesliceL1.rbegin(), back_inserter(*timesliceSN.rbegin())); } } if (parameters.writeL0()) { outputFile.put(*timeslice); } if (parameters.writeL1()) { outputFile.put(JTimesliceL1(timesliceL1, timesliceRouter, moduleRouter, parameters.TMaxLocal_ns)); } if (parameters.writeL2()) { outputFile.put(JTimesliceL1(timesliceL2, timesliceRouter, moduleRouter, parameters.L2.TMaxLocal_ns)); } if (parameters.writeSN()) { outputFile.put(JTimesliceL1(timesliceSN, timesliceRouter, moduleRouter, parameters.SN.TMaxLocal_ns)); } } STATUS(endl); JSingleFileScanner::typelist> io(inputFile); io >> outputFile; outputFile.close(); }