#include #include #include #include #include #include "km3net-dataformat/offline/Head.hh" #include "km3net-dataformat/offline/MultiHead.hh" #include "km3net-dataformat/offline/Evt.hh" #include "km3net-dataformat/definitions/fitparameters.hh" #include "JDAQ/JDAQEventIO.hh" #include "JDAQ/JDAQTimesliceIO.hh" #include "JDAQ/JDAQSummarysliceIO.hh" #include "JDetector/JDetector.hh" #include "JDetector/JDetectorToolkit.hh" #include "JDetector/JModuleRouter.hh" #include "JTrigger/JHit.hh" #include "JTrigger/JFrame.hh" #include "JTrigger/JTimeslice.hh" #include "JTrigger/JHitL0.hh" #include "JTrigger/JBuildL0.hh" #include "JTrigger/JTriggerParameters.hh" #include "JSupport/JSingleFileScanner.hh" #include "JSupport/JParallelFileScanner.hh" #include "JSupport/JFileRecorder.hh" #include "JSupport/JSupport.hh" #include "JSupport/JMeta.hh" #include "JFit/JPMTW0.hh" #include "JFit/JLine3Z.hh" #include "JFit/JGandalf.hh" #include "JFit/JLine3ZRegressor.hh" #include "JFit/JFitToolkit.hh" #include "JFit/JModel.hh" #include "JReconstruction/JHitW0.hh" #include "JReconstruction/JEvt.hh" #include "JReconstruction/JEvtToolkit.hh" #include "JPhysics/JConstants.hh" #include "JTools/JFunction1D_t.hh" #include "JTools/JFunctionalMap_t.hh" #include "JPhysics/JPDFTable.hh" #include "JPhysics/JPDFToolkit.hh" #include "Jeep/JParser.hh" #include "Jeep/JMessage.hh" /** * \file * * Program to perform JFIT::JRegressor fit with I/O of JFIT::JEvt data. * \author mdejong */ int main(int argc, char **argv) { using namespace std; using namespace JPP; using namespace KM3NETDAQ; typedef JParallelFileScanner< JTypeList > JParallelFileScanner_t; typedef JParallelFileScanner_t::multi_pointer_type multi_pointer_type; typedef JTYPELIST::typelist typelist; JParallelFileScanner_t inputFile; JFileRecorder outputFile; JLimit_t& numberOfEvents = inputFile.getLimit(); string detectorFile; string pdfFile; double roadWidth_m; double R_Hz; size_t numberOfPrefits; double TTS_ns; double E_GeV; int debug; try { JParser<> zap("Program to perform fit of muon path to data."); zap['f'] = make_field(inputFile); zap['o'] = make_field(outputFile) = "path.root"; zap['a'] = make_field(detectorFile); zap['n'] = make_field(numberOfEvents) = JLimit::max(); zap['P'] = make_field(pdfFile); zap['R'] = make_field(roadWidth_m) = numeric_limits::max(); zap['B'] = make_field(R_Hz) = 6.0e3; zap['N'] = make_field(numberOfPrefits) = 1; zap['T'] = make_field(TTS_ns); zap['E'] = make_field(E_GeV) = 1.0e3; 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 router(detector); typedef vector JDataL0_t; const JBuildL0 buildL0; typedef JRegressor JRegressor_t; JRegressor_t::debug = debug; JRegressor_t::T_ns.setRange(-50.0, +50.0); // ns JRegressor_t::Vmax_npe = 10.0; // npe JRegressor_t::MAXIMUM_ITERATIONS = 10000; JRegressor_t fit(pdfFile, TTS_ns); for (int i = 0; i != JRegressor_t::NUMBER_OF_PDFS; ++i) { fit.npe[i].transform(JRegressor_t::JNPE_t::transformer_type::getDefaultTransformer()); // get rid of R-dependent weight function } fit.estimator.reset(new JMEstimatorLorentzian()); fit.parameters.resize(4); fit.parameters[0] = JLine3Z::pX(); fit.parameters[1] = JLine3Z::pY(); fit.parameters[2] = JLine3Z::pDX(); fit.parameters[3] = JLine3Z::pDY(); if (fit.getRmax() < roadWidth_m) { roadWidth_m = fit.getRmax(); WARNING("Set road width to [m] " << roadWidth_m << endl); } const double Rmax_m = 100.0; // hit selection [m] const double Tmax_ns = 10.0; // hit selection [ns] outputFile.open(); outputFile.put(JMeta(argc, argv)); while (inputFile.hasNext()) { STATUS("event: " << setw(10) << inputFile.getCounter() << '\r'); DEBUG(endl); multi_pointer_type ps = inputFile.next(); const JDAQEvent* tev = ps; const JEvt* in = ps; JEvt cp = *in; JEvt out; cp.select(numberOfPrefits, qualitySorter); JDataL0_t dataL0; buildL0(*tev, router, true, back_inserter(dataL0)); for (JEvt::const_iterator track = cp.begin(); track != cp.end(); ++track) { const JRotation3D R (getDirection(*track)); const JLine1Z tz(getPosition (*track).rotate(R), track->getT()); const JModel match(tz, roadWidth_m, JRegressor_t::T_ns); // hit selection based on start value multiset top; JRange Z = JRange::DEFAULT_RANGE; for (JDataL0_t::const_iterator i = dataL0.begin(); i != dataL0.end(); ++i) { JHitW0 hit(*i, R_Hz); hit.rotate(R); if (match(hit)) { top.insert(hit.getPMTIdentifier()); const double t1 = hit.getT() - tz.getT(hit); if (tz.getDistance(hit) <= Rmax_m && t1 >= -Tmax_ns && t1 <= +Tmax_ns) { Z.include(hit.getZ() - tz.getDistance(hit) / getTanThetaC()); } } } vector buffer; for (JDetector::const_iterator module = detector.begin(); module != detector.end(); ++module) { JPosition3D pos(module->getPosition()); pos.rotate(R); if (tz.getDistance(pos) <= roadWidth_m && Z(pos.getZ() - tz.getDistance(pos) / getTanThetaC())) { for (unsigned int i = 0; i != module->size(); ++i) { const JDAQPMTIdentifier id(module->getID(), i); JPMT pmt(module->getPMT(i)); pmt.rotate(R); buffer.push_back(JPMTW0(pmt, R_Hz, top.count(id))); } } } const int NDF = buffer.size() - fit.parameters.size(); if (NDF >= 0) { // set fit parameters if (track->getE() > 0.1) fit.E_GeV = track->getE(); else fit.E_GeV = E_GeV; const double chi2 = fit(JLine3Z(tz), buffer.begin(), buffer.end()); JTrack3D tb(fit.value); tb.rotate_back(R); out.push_back(getFit(JMUONPATH, tb, getQuality(chi2), NDF)); out.rbegin()->setW(track->getW()); } } // apply default sorter sort(out.begin(), out.end(), qualitySorter); outputFile.put(out); } STATUS(endl); JSingleFileScanner::typelist> io(inputFile); io >> outputFile; outputFile.close(); }