#include #include #include #include #include #include "TROOT.h" #include "TFile.h" #include "TH1D.h" #include "TF1.h" #include "km3net-dataformat/offline/Head.hh" #include "km3net-dataformat/offline/Evt.hh" #include "km3net-dataformat/offline/Hit.hh" #include "km3net-dataformat/tools/time_converter.hh" #include "JDAQ/JDAQEventIO.hh" #include "JDAQ/JDAQTimesliceIO.hh" #include "JAAnet/JHead.hh" #include "JAAnet/JHeadToolkit.hh" #include "JAAnet/JAAnetToolkit.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/JSuperFrame2D.hh" #include "JTrigger/JHitL0.hh" #include "JTrigger/JHitL1.hh" #include "JTrigger/JBuildL0.hh" #include "JTrigger/JBuildL1.hh" #include "JTrigger/JBuildL2.hh" #include "JTrigger/JAlgorithm.hh" #include "JSupport/JTriggeredFileScanner.hh" #include "JSupport/JMonteCarloFileSupportkit.hh" #include "JSupport/JSupport.hh" #include "Jeep/JPrint.hh" #include "Jeep/JParser.hh" #include "Jeep/JMessage.hh" /** * \file * Example program to determine time slewing of L1 hits. * \author mdejong */ int main(int argc, char **argv) { using namespace std; using namespace JPP; using namespace KM3NETDAQ; JTriggeredFileScanner<> inputFile; JLimit_t& numberOfEvents = inputFile.getLimit(); string outputFile; string detectorFile; double Tmax_ns; double ctMin; int debug; try { JParser<> zap("Example program to determine time slewing of L1 hits."); zap['f'] = make_field(inputFile); zap['o'] = make_field(outputFile) = "hitL1.root"; zap['a'] = make_field(detectorFile); zap['n'] = make_field(numberOfEvents) = JLimit::max(); zap['T'] = make_field(Tmax_ns) = 15.0; // [ns] zap['C'] = make_field(ctMin) = 0.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 router(detector); detector -= get(getHeader(inputFile)); TFile out(outputFile.c_str(), "recreate"); vector H1D; const int nx = 100; const double xmin = -25.0; // [ns] const double xmax = +25.0; // [ns] for (int i = 0; i <= NUMBER_OF_PMTS; ++i) { ostringstream os; os << "M[" << setfill('0') << setw(2) << i << "]"; H1D.push_back(new TH1D(os.str().c_str(), NULL, nx, xmin, xmax)); } typedef vector JDataL0_t; typedef vector JDataL2_t; JSuperFrame2D buffer; // I/O buffer for time-sorted calibrated hits const JBuildL0 buildL0; const JBuildL2 buildL2(JL2Parameters(2, Tmax_ns, ctMin)); while (inputFile.hasNext()) { STATUS("event: " << setw(10) << inputFile.getCounter() << '\r'); DEBUG(endl); JTriggeredFileScanner<>::multi_pointer_type ps = inputFile.next(); const JDAQEvent* tev = ps; const Evt* event = ps; const time_converter converter(*event, *tev); vector::const_iterator muon = find_if(event->mc_trks.begin(), event->mc_trks.end(), is_muon); if (muon != event->mc_trks.end()) { // muon const JTrack3D trk = getTrack(*muon); JDataL0_t dataL0; JDataL2_t dataL2; buildL0(*tev, router, true, back_inserter(dataL0)); buildL2(*tev, router, true, back_inserter(dataL2)); // select first hit sort(dataL2.begin(), dataL2.end(), timeSorter); dataL2.erase(unique(dataL2.begin(), dataL2.end(), equal_to()), dataL2.end()); for (JDataL0_t::const_iterator hit = dataL0.begin(); hit != dataL0.end(); ++hit) { const double t0 = trk.getT(*hit); const double t1 = converter.getTime(hit->getT()); H1D[1]->Fill(t1 - t0); } for (JDataL2_t::const_iterator hit = dataL2.begin(); hit != dataL2.end(); ++hit) { const double t0 = trk.getT(*hit); const double t1 = converter.getTime(hit->begin()->getT()); const size_t index = min(hit->size(), H1D.size() - 1); H1D[index]->Fill(t1 - t0); } } } STATUS(endl); TF1 f1("f1", "[0]*exp(-0.5*(x-[1])*(x-[1])/([2]*[2])) + [3]"); string option = "LL"; if (debug < JEEP::debug_t && option.find('Q') == string::npos) { option += "Q"; } for (vector::iterator h1 = H1D.begin(); h1 != H1D.end(); ++h1) { // start values Double_t ymax = 0.0; Double_t x0 = 0.0; // [ns] Double_t sigma = 2.0; // [ns] Double_t ymin = 0.0; for (int i = 1; i != (*h1)->GetNbinsX(); ++i) { const Double_t x = (*h1)->GetBinCenter (i); const Double_t y = (*h1)->GetBinContent(i); if (y > ymax) { ymax = y; x0 = x; } } f1.SetParameter(0, ymax); f1.SetParameter(1, x0); f1.SetParameter(2, sigma); f1.FixParameter(3, ymin); // fit (*h1)->Fit(&f1, option.c_str(), "same", x0 - 3 * sigma, x0 + 3 * sigma); cout << "\tt0.push_back(" << showpos << FIXED(4,2) << f1.GetParameter(1) << ");" << endl; } out.Write(); out.Close(); }