#include #include #include #include "km3net-dataformat/offline/Head.hh" #include "km3net-dataformat/offline/MultiHead.hh" #include "km3net-dataformat/offline/Evt.hh" #include "km3net-dataformat/offline/Hit.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 "JDAQ/JDAQEvaluator.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/JHitL0.hh" #include "JTrigger/JHitL1.hh" #include "JTrigger/JBuildL1.hh" #include "JTrigger/JBuildL2.hh" #include "JTrigger/JTrigger3DShower.hh" #include "JTrigger/JTriggerMXShower.hh" #include "JTrigger/JTrigger3DMuon.hh" #include "JTrigger/JTriggerNB.hh" #include "JTrigger/JTriggerBits.hh" #include "JTrigger/JEventOverlap.hh" #include "JTrigger/JTimesliceRouter.hh" #include "JTrigger/JTriggeredEvent.hh" #include "JTrigger/JTimesliceL1.hh" #include "JTrigger/JTriggerParameters.hh" #include "JTrigger/JTriggerToolkit.hh" #include "JLang/JPipe.hh" #include "JROOT/JROOTClassSelector.hh" #include "JSupport/JSingleFileScanner.hh" #include "JSupport/JTreeScanner.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-trigger KM3NETDAQ::JDAQEvent 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; bool snapshot; JROOTClassSelection selection = getROOTClassSelection(); int debug; try { JParser<> zap("Auxiliary program to re-trigger event data."); zap['f'] = make_field(inputFile, "input file."); zap['o'] = make_field(outputFile, "output file.") = "trigger_reprocessor.root"; zap['n'] = make_field(numberOfEvents) = JLimit::max(); zap['a'] = make_field(detectorFile, "detector file."); zap['@'] = make_field(parameters, "trigger parameters") = JPARSER::initialised(); zap['U'] = make_field(reuse_parameters, "reuse trigger parameters from input file."); zap['S'] = make_field(snapshot, "use snapshot hits instead of triggered hits."); zap['C'] = make_field(selection, "selection of data types for output.") = JPARSER::initialised(); zap['d'] = make_field(debug, "debug flag.") = 1; zap(argc, argv); } catch(const exception &error) { FATAL(error.what() << endl); } setDAQLongprint(debug >= JEEP::debug_t); 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)); parameters.triggerNB.write.prescale = 1; DEBUG("Frame time [ms] " << getFrameTime() * 1e-6 << endl); DEBUG("Reset time [ms] " << getRTS() * 1e-6 << endl); DEBUG("Trigger" << endl << parameters << endl); const JModuleRouter moduleRouter(detector); if (parameters.writeSummary()) { WARNING("Discard writeSummary option during reprocesing of data." << endl); } if (parameters.writeL1()) { WARNING("Discard writeL1 option during reprocesing of data." << endl); } if (parameters.writeL2()) { WARNING("Discard writeL2 option during reprocesing of data." << endl); } if (parameters.writeSN()) { WARNING("Discard writeSN 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); const JBuildL2_t buildNB(parameters.NB); JTimesliceRouter timesliceRouter(parameters.numberOfBins); const JTriggerNB triggerNB (parameters); const JTrigger3DMuon trigger3DMuon (parameters); const JTrigger3DShower trigger3DShower(parameters); const JTriggerMXShower triggerMXShower(parameters, detector); outputFile.open(); if (!outputFile.is_open()) { FATAL("Error opening file " << outputFile << endl); } outputFile.put(JMeta(argc, argv)); outputFile.put(parameters); JTreeScanner::debug = debug; JTreeScanner scan(inputFile); if (scan.getEntries() == 0) { FATAL("No summary data." << endl); } while (inputFile.hasNext()) { STATUS("event: " << setw(10) << inputFile.getCounter() << '\r'); DEBUG(endl); JDAQEvent* evt = inputFile.next(); Long64_t index = scan.find(*evt); JDAQSummaryslice* summary = scan.getEntry(index); if (evt->getFrameIndex() != summary->getFrameIndex()) { DEBUG("Frame indices " << evt->getFrameIndex() << " != " << summary->getFrameIndex() << endl); } JDAQTimeslice timeslice(*evt, *summary, snapshot); 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()); JTimeslice_t timesliceNB(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)); // Nano-beacon trigger if (parameters.triggerNB.enabled) { JSuperFrame2D_t::iterator __end = partition(buffer.begin(), buffer.end(), parameters.triggerNB.pmts); if (buffer.begin() != __end) { timesliceNB.push_back(JSuperFrame1D_t(super_frame->getDAQChronometer(), super_frame->getModuleIdentifier(), module.getPosition())); JSuperFrame1D_t zbuf; buildL1(buffer.begin(), __end , back_inserter(zbuf)); buildNB(buffer.begin() , __end, zbuf, back_inserter(*timesliceNB.rbegin())); } } // 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())); } } // Trigger if (parameters.triggerNB.enabled) { const JTriggerInput trigger_input(timesliceNB); for (JTriggerInput::const_iterator hit = trigger_input.begin(); hit != trigger_input.end(); ++hit) { if (parameters.triggerNB.write()) { JTriggeredEvent tev(timesliceNB.getDAQChronometer(), getTriggerMask(triggerNB.getTriggerBit()), *hit, timesliceRouter, moduleRouter, parameters.TMaxLocal_ns, parameters.triggerNB.DMax_m, getTimeRange(parameters.triggerNB)); outputFile.put(tev); } } } JTriggerInput trigger_input(timesliceL2); JTriggerOutput trigger_output; trigger3DMuon (trigger_input, back_inserter(trigger_output)); trigger3DShower(trigger_input, back_inserter(trigger_output)); triggerMXShower(trigger_input, timesliceL0, back_inserter(trigger_output)); trigger_output.merge(JEventOverlap(parameters.TMaxEvent_ns)); DEBUG("Number of triggers: " << trigger_output.size() << endl); for (JTriggerOutput::const_iterator event = trigger_output.begin(); event != trigger_output.end(); ++event) { JTriggeredEvent tev(*event, timesliceRouter, moduleRouter, parameters.TMaxLocal_ns, getTimeRange(parameters)); tev.setCounter(evt->getCounter()); outputFile.put(tev); } } STATUS(endl); JSingleFileScanner io(inputFile); selection.remove(); io | JValve(selection) | outputFile; outputFile.close(); }