#include #include #include #include #include #include #include #include "TROOT.h" #include "TFile.h" #include "TH1D.h" #include "TH2D.h" #include "km3net-dataformat/offline/Head.hh" #include "km3net-dataformat/offline/Evt.hh" #include "km3net-dataformat/online/JDAQ.hh" #include "km3net-dataformat/online/JDAQPMTIdentifier.hh" #include "km3net-dataformat/online/JDAQSummaryFrame.hh" #include "JDAQ/JDAQEventIO.hh" #include "JDAQ/JDAQSummarysliceIO.hh" #include "JDetector/JDetector.hh" #include "JDetector/JDetectorToolkit.hh" #include "JDetector/JModuleRouter.hh" #include "JDetector/JPMTRouter.hh" #include "JTrigger/JHitL0.hh" #include "JTrigger/JBuildL0.hh" #include "JTrigger/JTriggerToolkit.hh" #include "JSupport/JSupport.hh" #include "JSupport/JTriggeredFileScanner.hh" #include "JSupport/JSummaryFileRouter.hh" #include "JSupport/JMonteCarloFileSupportkit.hh" #include "JAAnet/JAAnetToolkit.hh" #include "JFit/JLine1Z.hh" #include "JFit/JNPEHit.hh" #include "JFit/JModel.hh" #include "JPhysics/KM3NeT.hh" #include "JPhysics/JGeane.hh" #include "JReconstruction/JHitW0.hh" #include "JReconstruction/JEvt.hh" #include "JReconstruction/JEvtToolkit.hh" #include "JReconstruction/JMuonStartParameters_t.hh" #include "JReconstruction/JStart.hh" #include "JROOT/JManager.hh" #include "JTools/JAbstractHistogram.hh" #include "JLang/JPredicate.hh" #include "JLang/JComparator.hh" #include "Jeep/JPrint.hh" #include "Jeep/JParser.hh" #include "Jeep/JMessage.hh" namespace { using namespace JPP; using namespace KM3NETDAQ; /** * Get probability of given response in optical module due to random background. * * \param module module * \param frame summary frame * \param R_Hz multiples rates [Hz] * \param T_ns time window [ns] * \param top list of indentifiers of PMTs with hit in given module * \return probability */ inline double getProbability(const JModule& module, const JDAQSummaryFrame& frame, const JRateL1_t& R_Hz, const double T_ns, const std::multiset& top) { using namespace std; using namespace JPP; using namespace KM3NETDAQ; size_t N = 0; // number of active PMTs size_t M = 0; // multiplicity double R = 0.0; // total rate for (size_t i = 0; i != module.size(); ++i) { if (getDAQStatus(frame, module, i) && getPMTStatus(frame, module, i) && !module.getPMT(i).has(PMT_DISABLE)) { N += 1; M += top.count(i); R += frame.getRate(i); } } if (N > 0) return JFIT::getProbability(N, M, JK40Rates(R/N, R_Hz), T_ns); else return (M == 0 ? 1.0 : 0.0); } } /** * \file * * Program to test JMuonStart.cc. * * \author mdejong */ int main(int argc, char **argv) { using namespace std; using namespace JPP; using namespace KM3NETDAQ; typedef JTriggeredFileScanner JTriggeredFileScanner_t; typedef JTriggeredFileScanner_t::multi_pointer_type multi_pointer_type; typedef JAbstractHistogram JHistogram_t; JTriggeredFileScanner_t inputFile; JLimit_t& numberOfEvents = inputFile.getLimit(); string detectorFile; string outputFile; JMuonStartParameters_t parameters; JK40Rates rates_Hz; JHistogram_t histogram; int debug; try { parameters.numberOfPrefits = 1; JParser<> zap("Program to test JMuonStart."); zap['f'] = make_field(inputFile, "input file (output of JXXXMuonReconstruction.sh)"); zap['a'] = make_field(detectorFile); zap['n'] = make_field(numberOfEvents) = JLimit::max(); zap['o'] = make_field(outputFile) = "start.root"; zap['@'] = make_field(parameters) = JPARSER::initialised(); zap['B'] = make_field(rates_Hz) = KM3NET::getK40Rates(); zap['H'] = make_field(histogram, "histogram binning") = JHistogram_t(200, -250.0, +250.0); zap['d'] = make_field(debug) = 1; zap(argc, argv); } catch(const exception& error) { FATAL(error.what() << endl); } JDetector detector; try { load(detectorFile, detector); } catch(const JException& error) { FATAL(error); } const JModuleRouter moduleRouter(detector); const JPMTRouter pmtRouter (detector); JSummaryFileRouter summary(inputFile, parameters.R_Hz); const JTimeRange T_ns (parameters.TMin_ns, parameters.TMax_ns); const JStart start(parameters.Pmin1, parameters.Pmin2, parameters.Nmax2); typedef JHitL0 hit_type; typedef vector buffer_type; const JBuildL0 buildL0; JHead head; try { head = getHeader(inputFile); } catch(const JException& error) { FATAL(error); } const JPosition3D center = get(head); NOTICE("center: " << center << endl); JManager H1(new TH1D("[%].1L", NULL, histogram.getNumberOfBins(), histogram.getLowerLimit(), histogram.getUpperLimit())); JManager H2(new TH2D("[%].2D", NULL, histogram.getNumberOfBins(), histogram.getLowerLimit(), histogram.getUpperLimit(), histogram.getNumberOfBins(), histogram.getLowerLimit(), histogram.getUpperLimit())); while (inputFile.hasNext()) { STATUS("event: " << setw(10) << inputFile.getCounter() << '\r'); DEBUG(endl); multi_pointer_type ps = inputFile.next(); JDAQEvent* tev = ps; JEvt* evt = ps; Evt* event = ps; summary.update(*tev); vector::iterator muon = event->mc_trks.end(); for (vector::iterator track = event->mc_trks.begin(); track != event->mc_trks.end(); ++track) { if (is_muon(*track)) { if (muon == event->mc_trks.end() || track->E > muon->E) { muon = track; } } } if (muon == event->mc_trks.end()) { continue; } if (muon->len == 0.0) { muon->len = gWater(muon->E); } typedef JRange JRange_t; JRange_t za(JRange_t::DEFAULT_RANGE); // Monte Carlo JRange_t zb(JRange_t::DEFAULT_RANGE); // original result JRange_t zd(JRange_t::DEFAULT_RANGE); // new { const JRotation3D R(getDirection(*muon)); JTrack3D ta = getTrack(*muon); ta.add(center); ta.rotate(R); za.setLowerLimit(ta.getZ()); za.setLength(fabs(muon->len)); } buffer_type dataL0; buildL0(*tev, moduleRouter, true, back_inserter(dataL0)); if (!evt->empty()) { sort(evt->begin(), evt->end(), qualitySorter); JEvt::const_iterator track = evt->begin(); const JRotation3D R (getDirection(*track)); const JLine1Z tz(getPosition (*track).rotate(R), track->getT()); const JFIT::JModel match(tz, parameters.roadWidth_m, T_ns); if (track->has(JMUONSTART)) { if (track->hasW(JSTART_LENGTH_METRES) && track->getW(JSTART_LENGTH_METRES) > 0.0) { zb = JRange_t(tz.getZ(), tz.getZ() + track->getW(JSTART_LENGTH_METRES)); } } { map > top; for (buffer_type::const_iterator i = dataL0.begin(); i != dataL0.end(); ++i) { hit_type hit(*i); hit.rotate(R); if (match(hit)) { top[hit.getModuleID()].insert(hit.getPMTAddress()); } } struct JHit_t { double getZ() const { return z; } double getP() const { return p; } double z; double p; }; vector data; for (JDetector::const_iterator module = detector.begin(); module != detector.end(); ++module) { if (summary.hasSummaryFrame(module->getID())) { const JDAQSummaryFrame& frame = summary.getSummaryFrame(module->getID()); JPosition3D pos(module->getPosition()); pos.transform(R, tz.getPosition()); if (pos.getX() <= parameters.roadWidth_m) { const double z = pos.getZ() - pos.getX() / getTanThetaC(); const double p = getProbability(*module, frame, rates_Hz.getMultiplesRates(), T_ns.getLength(), top[module->getID()]); data.push_back({ z, p }); } } } if (!data.empty()) { sort(data.begin(), data.end(), make_comparator(&JHit_t::getZ)); vector::const_iterator track_start = start.find(data. begin(), data. end()); vector::const_reverse_iterator track_end = start.find(data.rbegin(), data.rend()); if (track_start == data. end()) { track_start = data. begin(); } if (track_end == data.rend()) { track_end = data.rbegin(); } zd = JRange_t(track_start->getZ(), track_end->getZ()); zd.add(tz.getZ()); } } } H1["B"]->Fill(zb.getLength() - za.getLength()); H1["D"]->Fill(zd.getLength() - za.getLength()); H2["B"]->Fill(zb.getLowerLimit() - za.getLowerLimit(), zb.getUpperLimit() - za.getUpperLimit()); H2["D"]->Fill(zd.getLowerLimit() - za.getLowerLimit(), zd.getUpperLimit() - za.getUpperLimit()); } STATUS(endl); TFile out(outputFile.c_str(), "recreate"); out << H1 << H2; out.Write(); out.Close(); }