Example program to test generation of KM3NETDAQ::JDAQTimeslice with random data.
More...
#include <string>
#include <iostream>
#include <iomanip>
#include <vector>
#include "TROOT.h"
#include "TFile.h"
#include "TH1D.h"
#include "TH2D.h"
#include "TRandom3.h"
#include "km3net-dataformat/online/JDAQ.hh"
#include "km3net-dataformat/online/JDAQClock.hh"
#include "km3net-dataformat/online/JDAQTimeslice.hh"
#include "JTimeslice/JRandomTimeslice.hh"
#include "JROOT/JManager.hh"
#include "Jeep/JParser.hh"
#include "Jeep/JMessage.hh"
Go to the source code of this file.
|
int | main (int argc, char **argv) |
|
Example program to test generation of KM3NETDAQ::JDAQTimeslice with random data.
- Author
- mdejong
Definition in file examples/JTimeslice/JRandomTimesliceWriter.cc.
◆ main()
int main |
( |
int |
argc, |
|
|
char ** |
argv |
|
) |
| |
Definition at line 29 of file examples/JTimeslice/JRandomTimesliceWriter.cc.
36 Long64_t numberOfSlices;
44 JParser<> zap(
"Example program to test generation of time slices with random data.");
48 zap[
'B'] =
make_field(rate_Hz,
"background rate [Hz]");
49 zap[
'N'] =
make_field(recycling,
"number of recycles / time interval for sampling data") = make_pair(0, 0.0);
55 catch(
const exception &error) {
56 FATAL(error.what() << endl);
60 gRandom->SetSeed(seed);
66 for (Long64_t counter = 0; counter != numberOfSlices; ++counter) {
68 STATUS(
"slice: " << setw(10) << counter <<
'\r');
DEBUG(endl);
81 for (
double t1 = 0.0 *
getFrameTime() + gRandom->Exp(period_ns); t1 < 0.5 *
getFrameTime(); t1 += gRandom->Exp(period_ns)) {
82 buffer.push_back(
JDAQHit(PMT_L0[0], t1, 0));
87 for (
double t1 = 0.5 *
getFrameTime() + gRandom->Exp(period_ns); t1 < 1.0 *
getFrameTime(); t1 += gRandom->Exp(period_ns)) {
88 buffer.push_back(
JDAQHit(PMT_L0[1], t1, 0));
92 timeslice[0].
add(buffer.size(), buffer.data());
97 h0->Fill((Double_t) hit->getT());
100 for (
size_t i = 1; i <= recycling.first; ++i) {
102 timeslice.
recycle(recycling.second);
107 h0->Fill((Double_t) hit->getT());
#define DEBUG(A)
Message macros.
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Utility class to parse command line options.
Auxiliary class to manage set of compatible ROOT objects (e.g. histograms) using unique keys.
unsigned char JPMT_t
PMT channel in FPGA.
JDAQTimeslice & add(const JDAQTimeslice ×lice)
Add another timeslice.
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
KM3NeT DAQ data structures and auxiliaries.
double getFrameTime()
Get frame time duration.
static const int NUMBER_OF_PMTS
Total number of PMTs in module.
Timeslice with random data.
void recycle(const double T_ns)
Recycle time slice by randomly shuffling time intervals of data.