#include #include #include #include #include #include #include "TROOT.h" #include "TFile.h" #include "TH1D.h" #include "TF1.h" #include "TMath.h" #include "km3net-dataformat/online/JDAQ.hh" #include "JDAQ/JDAQEventIO.hh" #include "JDAQ/JDAQTimesliceIO.hh" #include "JDAQ/JDAQEvaluator.hh" #include "JDetector/JDetector.hh" #include "JDetector/JDetectorToolkit.hh" #include "JDetector/JDAQHitRouter.hh" #include "JDetector/JPMTParametersMap.hh" #include "JTrigger/JBuildL0.hh" #include "JTrigger/JBuildL1.hh" #include "JLang/JSinglePointer.hh" #include "JROOT/JManager.hh" #include "JTools/JWeight.hh" #include "JROOT/JROOTClassSelector.hh" #include "JSupport/JSingleFileScanner.hh" #include "JSupport/JTreeScannerInterface.hh" #include "JSupport/JTreeScanner.hh" #include "JSupport/JAutoTreeScanner.hh" #include "JSupport/JSupport.hh" #include "JSupport/JMeta.hh" #include "Jeep/JPrint.hh" #include "Jeep/JParser.hh" #include "Jeep/JMessage.hh" namespace { // status of module enum JStatus_t { DEFAULT = 0, //!< module exists in detector NODATA = -2, //!< module not present in data or no entries OUT_SYNC = -3, //!< module is out-of-sync ERROR = -4, //!< module is out-of-sync, possibly multiple time offsets IN_SYNC = +1 //!< module is in-sync }; } /** * \file * * Example program to search for correlations between triggered events and timeslice data. * \author mdejong */ int main(int argc, char **argv) { using namespace std; using namespace JPP; using namespace KM3NETDAQ; JSingleFileScanner<> inputFile; JLimit_t& numberOfEvents = inputFile.getLimit(); string outputFile; string detectorFile; double TMaxLocal_ns; string pmtFile; int numberOfTimeslices; double binWidth_ns; double Pmin; JROOTClassSelector selector; int qaqc; int debug; string peaks; try { JParser<> zap("Example program to search for correlations between triggered events and timeslice data."); zap['f'] = make_field(inputFile); zap['o'] = make_field(outputFile) = ""; zap['a'] = make_field(detectorFile); zap['n'] = make_field(numberOfEvents) = JLimit::max(); zap['T'] = make_field(TMaxLocal_ns) = 20.0; zap['P'] = make_field(pmtFile) = ""; zap['N'] = make_field(numberOfTimeslices) = 100; zap['W'] = make_field(binWidth_ns) = 10e3; zap['p'] = make_field(Pmin) = 1.0e-6; zap['C'] = make_field(selector) = "", getROOTClassSelection(); zap['Q'] = make_field(qaqc) = 0; zap['d'] = make_field(debug) = 2; zap['O'] = make_field(peaks) = ""; //Make a small summary file, with the list of the peaks zap.read(argc, argv); } catch(const exception& error) { FATAL(error.what() << endl); } JDetector detector; try { load(detectorFile, detector); } catch(const JException& error) { FATAL(error); } const JDAQHitRouter router(detector); const double TMax_ns = max(binWidth_ns, getMaximalTime(detector)); typedef double hit_type; JBuildL0 buildL0; JBuildL1 buildL1(JBuildL1Parameters(TMaxLocal_ns, true)); vector data; typedef JManager JManager_t; const Double_t xmin = -(numberOfTimeslices + 1) * getFrameTime() - 0.5*binWidth_ns; const Double_t xmax = +(numberOfTimeslices + 1) * getFrameTime() + 0.5*binWidth_ns; const Int_t nx = (Int_t) ((xmax - xmin) / binWidth_ns); JManager_t manager(new TH1D("M_%", NULL, nx, xmin, xmax)); JTreeScanner<>::debug = debug; JSinglePointer< JTreeScannerInterface > ps; JAutoTreeScanner zmap = JType(); if (selector == "") { if ((ps = new JTreeScanner< JAssertConversion >(inputFile))->getEntries() != 0 || (ps = new JTreeScanner< JAssertConversion >(inputFile))->getEntries() != 0 || (ps = new JTreeScanner< JAssertConversion >(inputFile))->getEntries() != 0 || (ps = new JTreeScanner< JAssertConversion >(inputFile))->getEntries() != 0 || (ps = new JTreeScanner< JAssertConversion >(inputFile))->getEntries() != 0) { } else { FATAL("No timeslice data." << endl); } NOTICE("Selected class " << ps->getClassName() << endl); } else { ps = zmap[selector]; ps->configure(inputFile); } /* Updated method : * OOS DOMs can have an impact on the histogram of other DOMs * especialy OOS DOMs on adjacent In Sync or OOS DOMs (presence of ghost peaks) * To remove these ghost peaks, the idea is to iterate the method, removing * at each loop the most OOS DOM, and recomputing the histograms without * them (removing events where they are involved) * Once there is no more OOS DOMs in the remaining DOMs, end of precedure. * Steps of procedure (from J. Hofestaedt, code updated by R. Le Breton) * 1) Check if any DOM shows a peak away from 0 (OOS) * If no OOS -> finish procedure * 2) If at least one OOS, then identify the DOM that is most OOS, * ie. highest peak (other than delta time = 0) * 3) Remove the most OOS DOM and calculate the histograms again but do not * consider the OOS DOM (in later iterations more than one DOM is OOS and not * considered for histogram filling). * 4) start with 1) again */ ps->setLimit(inputFile.getLimit()); typedef map > map_type; // frame index -> event times bool OOS = true; //Boolean: to continue or not the procedure map_type buffer; map_type peaksPerDoms; map status2; //OOS DOMs IDs : if at least one OOS DOM, one ID is added into oosDomsId at each iteration map oosDomsId; //Default value for all modules : false for (JDetector::const_iterator module = detector.begin(); module != detector.end(); ++module) { oosDomsId[module->getID()] = false; } do{ for (JTreeScanner in(inputFile.getFilename()); in.hasNext(); ) { STATUS("event: " << setw(10) << in.getCounter() << '\r'); DEBUG(endl); JDAQEvent* event = in.next(); // Average time of triggered hits double t0 = 0.0; bool test_OOS = false; //If an OOS DOM is involved in an event, remove this event for (JDAQEvent::const_iterator hit = event->begin(); hit != event->end(); ++hit) { t0 += getTime(*hit, router.getPMT(*hit)); if(oosDomsId[hit->getModuleID()]) { test_OOS = true; break; } } if(test_OOS) continue; //If true, event not taken into account t0 /= event->size(); t0 += event->getFrameIndex() * getFrameTime(); buffer[event->getFrameIndex()].push_back(t0); } STATUS(endl); //Events : start from the beginning ps->rewind(); while (ps->hasNext()) { STATUS("event: " << setw(10) << ps->getCounter() << '\r'); DEBUG(endl); JDAQTimeslice* timeslice = ps->next(); map_type::const_iterator p = buffer.lower_bound(timeslice->getFrameIndex() - numberOfTimeslices); map_type::const_iterator q = buffer.upper_bound(timeslice->getFrameIndex() + numberOfTimeslices); int number_of_events = 0; for (map_type::const_iterator i = p; i != q; ++i) { number_of_events += i->second.size(); } if (number_of_events == 0) { continue; } for (JDAQTimeslice::const_iterator frame = timeslice->begin(); frame != timeslice->end(); ++frame) { data.clear(); buildL1(*frame, router.getModule(frame->getModuleID()), back_inserter(data)); TH1D* h1 = manager[frame->getModuleID()]; for (vector::const_iterator hit = data.begin(); hit != data.end(); ++hit) { const double t1 = *hit + frame->getFrameIndex() * getFrameTime(); for (map_type::const_iterator i = p; i != q; ++i) { for (map_type::mapped_type::const_iterator j = i->second.begin(); j != i->second.end(); ++j) { const double t0 = *j; h1->Fill(t1 - t0); } } } } } buffer.clear(); //Clear the buffer for the next loop, if any STATUS(endl); Double_t most_OOS_peak = 0; Int_t most_OOS_ID = 0; TF1 f1("f1", "[0]*exp(-0.5*(x-[1])*(x-[1])/([2]*[2])) + [3]"); string option = "L"; if (debug < debug_t) { option += "Q"; } for (JDetector::const_iterator module = detector.begin(); module != detector.end(); ++module) { bool status = false; status2[module->getID()] = DEFAULT; JManager_t::iterator p = manager.find(module->getID()); if (p == manager.end() || p->second->GetEntries() == 0) { status2[module->getID()] = NODATA; continue; } TH1D* h1 = p->second; // start values Double_t ymax = 0.0; Double_t x0 = 0.0; // [ns] Double_t sigma = 250.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; } ymin += y; } ymin /= h1->GetNbinsX(); f1.SetParameter(0, ymax); f1.SetParameter(1, x0); if (binWidth_ns < sigma) f1.SetParameter(2, sigma); else f1.FixParameter(2, binWidth_ns/sqrt(12.0)); f1.SetParameter(3, ymin); // fit h1->Fit(&f1, option.c_str(), "same", x0 - 5 * sigma, x0 + 5 * sigma); status = (fabs(f1.GetParameter(1)) <= 0.5*getFrameTime() && // position f1.GetParameter(0) >= f1.GetParameter(3)); // S/N if (status) status2[module->getID()] = IN_SYNC; // look for peaks at offsets equal to a multiple of frame time JWeight bg; // background map sn; // signal(s) for (Int_t i = 1, ns = -(numberOfTimeslices + 1); i <= h1->GetNbinsX(); ++i) { const Double_t x = h1->GetBinCenter (i); const Double_t y = h1->GetBinContent(i); while (x > (ns + 1) * getFrameTime() - TMax_ns) { ++ns; } if (fabs(x - ns * getFrameTime()) < TMax_ns) sn[ns].put(y); else bg .put(y); } DEBUG("Module " << setw(8) << module->getID() << endl); for (map::iterator i = sn.begin(); i != sn.end(); ) { const Double_t y = bg.getTotal() * i->second.getCount() / bg.getCount(); const Double_t P = TMath::PoissonI(i->second.getTotal(), y); DEBUG("Peak at " << setw(4) << i->first << " [frame time] " << "S/N = " << FIXED(7,1) << i->second.getTotal() << " / " << FIXED(7,1) << y << ' ' << SCIENTIFIC(12,3) << P << endl); if (P < Pmin && i->second.getTotal() > y) //Avoid edge effect (second bolean) ++i; // keep else sn.erase(i++); // remove } if (!(sn.size() == 1 && sn.begin()->first == 0)) { WARNING("Module " << setw(8) << module->getID() << " number of peaks " << sn.size() << endl); for (map::const_iterator i = sn.begin(); i != sn.end(); ++i) { const Double_t noise = bg.getTotal() * i->second.getCount() / bg.getCount(); WARNING("Peak at " << setw(4) << i->first << " [frame time] " << "S/N = " << FIXED(7,1) << i->second.getTotal() << " / " << FIXED(7,1) << noise << endl); peaksPerDoms[module->getID()].push_back(sn.size()); peaksPerDoms[module->getID()].push_back(i->first); peaksPerDoms[module->getID()].push_back(i->second.getTotal()); peaksPerDoms[module->getID()].push_back(noise); } status2[module->getID()] = (sn.size() == 1 ? OUT_SYNC : ERROR); } //if DOM is OOS if (!status) { //if OOS peak greater than the last and not in OOS DOMs vector if ((f1.GetParameter(0) > most_OOS_peak) && !(oosDomsId[module->getID()])){ most_OOS_peak = f1.GetParameter(0); most_OOS_ID = module->getID(); } } } if (most_OOS_ID != 0) { oosDomsId[most_OOS_ID] = true;//exclude the DOM to compute the histograms //Reseting all histograms for (JDetector::const_iterator module = detector.begin(); module != detector.end(); ++module) { manager[module->getID()]->Reset("ICESM"); } peaksPerDoms.clear(); } else{ OOS = false; } }while(OOS); //To manage non working and OOS DOMs NOTICE("Managing non-working and OOS DOMs." << endl); if (pmtFile != "") { JPMTParametersMap parameters; try { parameters.load(pmtFile.c_str()); } catch(const JException& error) {} for (map::const_iterator i = status2.begin(); i != status2.end(); ++i) { if (i->second != IN_SYNC) { NOTICE("Module " << setw(8) << i->first << " set QEs of all PMTs to zero." << endl); for (int pmt = 0; pmt != NUMBER_OF_PMTS; ++pmt) { parameters[JPMTIdentifier(i->first, pmt)].QE = 0.0; } } } ofstream out(pmtFile.c_str()); parameters.comment.add(JMeta(argc, argv)); out << parameters << endl; out.close(); } if (outputFile != "") { manager.Write(outputFile.c_str()); } //print small summary file of the peaks if (peaks != ""){ ofstream stream(peaks); stream << "#DOM NB_PEAKS TIMESLICE_SHIFT SIGNAL NOISE\n"; for (map_type::const_iterator dom = peaksPerDoms.begin(); dom != peaksPerDoms.end(); ++dom) { stream << dom->first << " "; int i = 1; for (map_type::mapped_type::const_iterator data = dom->second.begin(); data != dom->second.end(); ++data) { stream << *data; ((i%4 == 0) ? stream << "\n" : stream << " "); i++; } } stream.close(); } //QAQC int nin = 0; int nout = 0; for (map::const_iterator i = status2.begin(); i != status2.end(); ++i) { if (i->second == IN_SYNC) { ++nin; } if (i->second != IN_SYNC && i->second != NODATA) { ++nout; } } NOTICE("Number of modules out/in sync " << nout << '/' << nin << endl); QAQC(nin << ' ' << nout << endl); return 0; }