#include #include #include #include #include #include "JDAQ/JDAQ.hh" #include "JDAQ/JDAQClock.hh" #include "JDAQ/JDAQTimeslice.hh" #include "JDAQ/JDAQEvent.hh" #include "JDAQ/JDAQSummaryslice.hh" #include "JDetector/JDetector.hh" #include "JDetector/JDetectorToolkit.hh" #include "JDetector/JModuleRouter.hh" #include "JTrigger/JHit.hh" #include "JTrigger/JHitToolkit.hh" #include "JTrigger/JFrame.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/JTrigger3DMuon.hh" #include "JTrigger/JEventOverlap.hh" #include "JTrigger/JTimesliceRouter.hh" #include "JTrigger/JTriggeredEvent.hh" #include "JTrigger/JTimesliceL1.hh" #include "JTrigger/JTriggerParameters.hh" #include "JLang/JConversionIterator.hh" #include "JSupport/JMultipleFileScanner.hh" #include "JSupport/JFileRecorder.hh" #include "JSupport/JSupport.hh" #include "Jeep/JParser.hh" #include "Jeep/JMessage.hh" #include "Jeep/JTimer.hh" /** * Auxiliary program to trigger timeslice data. */ int main(int argc, char **argv) { using namespace std; using namespace KM3NETDAQ; using namespace JSUPPORT; using namespace JTRIGGER; JMultipleFileScanner<> inputFile; JFileRecorder outputFile; JLimit_t& numberOfEvents = inputFile.getLimit(); string detectorFile; JTriggerParameters parameters; set triggerBits; bool reprocess; JDAQClock clock; int debug; try { JParser<> zap; zap['f'] = make_field(inputFile); zap['o'] = make_field(outputFile) = "trigger_processor.dat"; zap['n'] = make_field(numberOfEvents) = JLimit::max(); zap['a'] = make_field(detectorFile); zap['@'] = make_field(parameters); zap['T'] = make_field(triggerBits); zap['R'] = make_field(reprocess); zap['c'] = make_field(clock) = JDAQClock::KM3NET, JDAQClock::PPM_DU, JDAQClock::ANTARES; zap['d'] = make_field(debug) = 0; zap['@'] = JPARSER::initialised(); zap['T'] = JPARSER::initialised(); if (zap.read(argc, argv) != 0) return 1; } catch(const exception &error) { FATAL(error.what() << endl); } using namespace JSIRENE; using namespace KM3NETDAQ; using JEEP::JTimer; using JTOOLS::make_range; using JLANG::JConversionIterator; cout.tie(&cerr); setLongprint(cout, debug >= 3); DEBUG("Frame time [ms] " << getFrameTime() * 1e-6 << endl); DEBUG("Reset time [ms] " << getRTS() * 1e-6 << endl); JDetector detector; try { load(detectorFile, detector); } catch(const JException& error) { FATAL(error); } parameters.set(getMaximalDistance(detector)); const JTrigger3DMuon trigger3DMuon (parameters); const JTrigger3DShower trigger3DShower(parameters); const JModuleRouter moduleRouter(detector); //typedef JHit hit_t; //typedef int hit_t; typedef double hit_t; typedef vector JFrameL1_t; typedef JSuperFrame JSuperFrame_t; typedef JTimeslice JTimesliceL1_t; typedef JBuildL1 JBuildL1_t; typedef JBuildL2 JBuildL2_t; const JBuildL1_t buildL1(JHitToolkit::getT(parameters.TMaxLocal_ns)); JTimer timerCC("Calibration"); JTimer timerL1("L1"); JTimer timerL2("L2"); JTimer timerRX("Timeslice router"); JTimer timerTR("Trigger"); JTimer timerTX("Trigger router"); JTimer timerTW("Timeslice writer"); JTimer timerSW("Summary writer"); outputFile.open(); if (!outputFile.is_open()) FATAL("Error opening file " << outputFile << endl); outputFile.put(parameters); unsigned int numberOfTriggers = 0; JObjectIterator* input = NULL; if (!reprocess) input = new JMultipleFileScanner(inputFile); else input = new JConversionIterator(new JMultipleFileScanner(inputFile)); for (JCounter_t counter = 0; input->hasNext(); ++counter) { STATUS("event: " << setw(10) << counter << '\r'); DEBUG(endl); JDAQTimeslice* timeslice = input->next(); DEBUG(*timeslice << endl); timerRX.start(); const JTimesliceRouter timesliceRouter(*timeslice, parameters.numberOfBins); timerRX.stop(); JFrameL1_t frameL1; JTimesliceL1_t timesliceL1(timeslice->getDAQChronometer()); JSuperFrame_t buffer; for (JDAQTimeslice::const_iterator super_frame = timeslice->begin(); super_frame != timeslice->end(); ++super_frame) { if (!moduleRouter.hasModule(super_frame->getModuleID())) { if (counter%100 == 0) ERROR("Missing module: " << super_frame->getModuleID() << endl); continue; } // calibration timerCC.start(); const JModule& module = moduleRouter.getModule(super_frame->getModuleID()); buffer(*super_frame, module); timerCC.stop(); // Apply high-rate veto buffer.applyHighRateVeto(parameters.highRateVeto_Hz); // L1 timerL1.start(); frameL1.clear(); buildL1(buffer, back_inserter(frameL1)); timerL1.stop(); // L2 timerL2.start(); const JBuildL2_t buildL2(buffer, parameters.L2Min, parameters.TMaxLocal_ns, parameters.ctMin); timesliceL1.push_back(JFrame(super_frame->getDAQChronometer(), JDAQPMTIdentifier(super_frame->getModuleID(), -1), module.getPosition())); buildL2(frameL1, back_inserter(timesliceL1.rbegin()->getData())); timerL2.stop(); } // Trigger timerTR.start(); JTriggerInput trigger_input(timesliceL1); JTriggerOutput trigger_output; if (triggerBits.count(trigger3DMuon .getTriggerBit())) trigger3DMuon (trigger_input, back_inserter(trigger_output)); if (triggerBits.count(trigger3DShower.getTriggerBit())) trigger3DShower(trigger_input, back_inserter(trigger_output)); trigger_output.merge(JEventOverlap(parameters.TMaxEvent_ns)); numberOfTriggers += trigger_output.size(); timerTR.stop(); DEBUG("Number of triggers: " << trigger_output.size() << endl); for (JTriggerOutput::const_iterator event = trigger_output.begin(); event != trigger_output.end(); ++event) { timerTX.start(); JTriggeredEvent tev(*event, timesliceRouter, moduleRouter, parameters.TMaxLocal_ns, parameters.getTimeRange()); timerTX.stop(); outputFile.put(tev); } if (parameters.writeL1) { timerTW.start(); outputFile.put(JTimesliceL1(timesliceL1, timesliceRouter, moduleRouter, parameters.TMaxLocal_ns)); timerTW.stop(); } if (parameters.writeSummary) { timerSW.start(); outputFile.put(JDAQSummaryslice(*timeslice)); timerSW.stop(); } } STATUS(endl); if (!reprocess && debug > 0) { const JCounter_t counter = dynamic_cast*>(input)->getCounter(); const double factor = 1.0 / (double) counter; for (const JTimer* buffer[] = { &timerCC, &timerL1, &timerL2, &timerRX, &timerTR, &timerTX, &timerTW, &timerSW, NULL }, **i = buffer; *i != NULL; ++i) (*i)->print(cout, factor); cout << "Number of trigger/slices " << numberOfTriggers << "/" << counter << endl; cout << "Trigger rate [Hz] " << numberOfTriggers * 1.0e9 * factor / getFrameTime() << endl; } if (input != NULL) delete input; outputFile.close(); }