#include #include #include #include #include #include #include "TROOT.h" #include "TFile.h" #include "TH1D.h" #include "km3net-dataformat/online/JDAQClock.hh" #include "JDAQ/JDAQTimesliceIO.hh" #include "JSupport/JSingleFileScanner.hh" #include "JSupport/JSupportToolkit.hh" #include "JSupport/JSupport.hh" #include "JROOT/JROOTClassSelector.hh" #include "JROOT/JManager.hh" #include "JLang/JObjectMultiplexer.hh" #include "JTrigger/JChecksum.hh" #include "Jeep/JParser.hh" #include "Jeep/JMessage.hh" namespace { /** * Print DAQ hit. * * \param out output stream * \param hit DAQ hit */ inline void print(std::ostream& out, const KM3NETDAQ::JDAQHit& hit) { using namespace std; out << "\thit: " << setw(3) << (int) hit.getPMT() << ' ' << setw(8) << hex << hit.getT() << ' ' << setw(10) << dec << hit.getT() << ' ' << setw(3) << (int) hit.getToT(); } /** * Print error type. * * \param out output stream * \param type error type */ inline void print(std::ostream& out, const int type) { using namespace std; using namespace JPP; switch (type) { case JChecksum::TIME_t: out << "time "; break; case JChecksum::ETDC_t: out << "TDC "; break; case JChecksum::EPMT_t: out << "PMT "; break; case JChecksum::EUDP_t: out << "UDP "; break; case JChecksum::SIZE_t: out << "size "; break; default: out << "? "; break; } } } /** * \file * * Example program to check KM3NETDAQ::JDAQTimeslice for errors. * \author mdejong */ int main(int argc, char **argv) { using namespace std; using namespace JPP; using namespace KM3NETDAQ; JSingleFileScanner inputFile; JLimit_t & numberOfEvents = inputFile.getLimit(); JROOTClassSelector selector; string outputFile; int N; double rate_Hz; int debug; try { JParser<> zap("Example program to check timeslice data for errors."); zap['f'] = make_field(inputFile); zap['n'] = make_field(numberOfEvents) = JLimit::max(); zap['C'] = make_field(selector) = getROOTClassSelection(); zap['o'] = make_field(outputFile) = "checksum.root"; zap['N'] = make_field(N, "number of rows to print") = 10; zap['R'] = make_field(rate_Hz, "high-rate veto [Hz]") = 0.0; zap['d'] = make_field(debug) = 3; zap(argc, argv); } catch(const exception& error) { FATAL(error.what() << endl); } JObjectMultiplexer in(inputFile, selector); if (rate_Hz > 0.0) { MAXIMAL_FRAME_SIZE = (int) (getFrameTime() * 1.0e-9 * rate_Hz * NUMBER_OF_PMTS + 0.5); } map errors; JFrameIndexRange frame_index = getFrameIndexRange(inputFile); const Long64_t nx = frame_index.second - frame_index.first + 1; const Double_t xmin = (Double_t) frame_index.first - 0.5; const Double_t xmax = (Double_t) frame_index.second + 0.5; JManager H1(new TH1D("[%]", NULL, nx, xmin, xmax)); for (counter_type counter = 0; in.hasNext(); ++counter) { STATUS("event: " << setw(10) << counter << '\r'); DEBUG(endl); const JDAQTimeslice* timeslice = in.next(); for (JDAQTimeslice::const_iterator frame = timeslice->begin(); frame != timeslice->end(); ++frame) { const JChecksum::result_type& result = checksum(*frame); for (JChecksum::const_iterator error = result.begin(); error != result.end(); ++error) { errors[error->type] += 1; H1[error->type]->Fill((Double_t) frame->getFrameIndex()); } if (debug >= debug_t) { for (JChecksum::const_iterator error = result.begin(); error != result.end(); ++error) { if (error->pos >= 0) { const int pmt = (*frame)[error->pos].getPMT(); cout << "Module " << setw(10) << frame->getModuleID() << '.' << setw(2) << setfill('0') << pmt << setfill(' ') << ' ' << setw(8) << frame->getFrameIndex() << ' ' << setw(6) << error->pos << '/' << setw(8) << frame->size() << ' ' << setw(2) << error->type << ' '; if (pmt < NUMBER_OF_PMTS) { cout << setw(1) << frame->testHighRateVeto(pmt) << setw(1) << frame->testFIFOStatus (pmt) << ' '; } cout << setw(4) << frame->getUDPNumberOfReceivedPackets() << '/' << setw(4) << frame->getUDPMaximalSequenceNumber(); cout << " "; print(cout, error->type); cout << endl << endl; deque buffer; for (int i = error->pos - 1, n = 0; i >= 0 && n <= N; --i) { if ((*frame)[i].getPMT() == (*frame)[error->pos].getPMT()) { buffer.push_front((*frame)[i]); ++n; } } for (deque::const_iterator i = buffer.begin(); i != buffer.end(); ++i) { print(cout, *i); cout << endl; } print(cout, (*frame)[error->pos]); cout << " <<< " << endl; for (int i = error->pos + 1, n = 0; i < frame->size() && n <= N; ++i) { if ((*frame)[i].getPMT() == (*frame)[error->pos].getPMT()) { print(cout, (*frame)[i]); cout << endl; ++n; } } } else { cout << "Module " << setw(10) << frame->getModuleID() << ' ' << setw(2) << ' ' << ' ' << setw(8) << frame->getFrameIndex() << ' ' << setw(6) << error->pos << '/' << setw(8) << frame->size() << ' ' << setw(2) << error->type << ' '; cout << setw(1) << ' ' << setw(1) << ' ' << ' '; cout << setw(4) << frame->getUDPNumberOfReceivedPackets() << '/' << setw(4) << frame->getUDPMaximalSequenceNumber(); cout << " "; print(cout, error->type); cout << endl << endl; } } } } } STATUS(endl); for (map::const_iterator i = errors.begin(); i != errors.end(); ++i) { print(cout, i->first); cout << ' ' << setw(8) << i->second << endl; } H1.Write(outputFile.c_str()); }