Jpp  debug
the software that should make you happy
Functions
software/JTimeslice/JRandomTimesliceWriter.cc File Reference

Auxiliary program to write KM3NETDAQ::JDAQTimeslice with random data. More...

#include <string>
#include <iostream>
#include <fstream>
#include "TRandom3.h"
#include "km3net-dataformat/offline/Head.hh"
#include "km3net-dataformat/online/JDAQClock.hh"
#include "JAAnet/JHead.hh"
#include "JAAnet/JHeadToolkit.hh"
#include "JDetector/JDetector.hh"
#include "JDetector/JDetectorToolkit.hh"
#include "JDetector/JPMTParametersMap.hh"
#include "JDetector/JK40DefaultSimulator.hh"
#include "JDetector/JPMTDefaultSimulator.hh"
#include "JDetector/JCLBDefaultSimulator.hh"
#include "JDetector/JDetectorSimulator.hh"
#include "JDAQ/JDAQSummarysliceIO.hh"
#include "JDAQ/JDAQTimesliceIO.hh"
#include "JTimeslice/JRandomTimeslice.hh"
#include "JTrigger/JSummaryRouter.hh"
#include "JTrigger/JK40RunByRunSimulator.hh"
#include "JTrigger/JPMTRunByRunSimulator.hh"
#include "JSupport/JSingleFileScanner.hh"
#include "JSupport/JFileRecorder.hh"
#include "JSupport/JSupport.hh"
#include "JSupport/JSupportToolkit.hh"
#include "JSupport/JRunByRun.hh"
#include "JSupport/JMeta.hh"
#include "Jeep/JTimer.hh"
#include "Jeep/JParser.hh"
#include "Jeep/JMessage.hh"

Go to the source code of this file.

Functions

int main (int argc, char **argv)
 

Detailed Description

Auxiliary program to write KM3NETDAQ::JDAQTimeslice with random data.


To pipe data to application JTriggerProcessor.cc whilst maintaining header and meta data, two output files can be specified.
The file names should have extensions ".root" and ".dat", respectively.
The first corresponds to a ROOT formatted file with header and meta data and the second to the pipe with timeslice data.

Author
mdejong

Definition in file software/JTimeslice/JRandomTimesliceWriter.cc.

Function Documentation

◆ main()

int main ( int  argc,
char **  argv 
)

Definition at line 48 of file software/JTimeslice/JRandomTimesliceWriter.cc.

49 {
50  using namespace std;
51  using namespace JPP;
52  using namespace KM3NETDAQ;
53 
55 
57  string detectorFile;
58  Long64_t numberOfSlices;
59  JDAQHit::JTDC_t TCLB_ns;
60  JPMTParametersMap pmtParameters;
61  JK40Rates rates_Hz;
62  int run;
63  JRunByRun runbyrun;
64  pair<size_t, double> recycling;
65  double sigma_ns;
66  bool fast;
67  UInt_t seed;
68  int debug;
69 
70  try {
71 
72  JParser<> zap("Auxiliary program to write time slices with random data.");
73 
74  zap['o'] = make_field(outputFile, "output file");
75  zap['n'] = make_field(numberOfSlices);
76  zap['a'] = make_field(detectorFile, "detector.");
77  zap['R'] = make_field(run, "run number") = -1;
78  zap['r'] = make_field(runbyrun, "option for run-by-run mode") = JPARSER::initialised();
79  zap['P'] = make_field(pmtParameters, "PMT simulation data (or corresponding file name)") = JPARSER::initialised();
80  zap['B'] = make_field(rates_Hz, "background rates [Hz]") = JPARSER::initialised();
81  zap['T'] = make_field(TCLB_ns, "CLB state-machine time jitter") = 256; // [ns]
82  zap['N'] = make_field(recycling, "number of recycles / time interval for sampling data [ns]") = make_pair(0, 0.0);
83  zap['s'] = make_field(sigma_ns, "intrinsic time smearing of K40 coincidences [ns]") = JK40DefaultSimulatorInterface::getSigma();
84  zap['F'] = make_field(fast, "fast - disable PMT simulation");
85  zap['S'] = make_field(seed, "seed") = 0;
86  zap['d'] = make_field(debug, "debug") = 0;
87 
88  zap(argc, argv);
89  }
90  catch(const exception &error) {
91  FATAL(error.what() << endl);
92  }
93 
94 
95  if (outputFile.size() != 1 &&
96  outputFile.size() != 2) {
97  FATAL("Invalid number of output files; should be 1 or 2." << endl);
98  }
99 
100  gRandom->SetSeed(seed);
101 
102  JK40DefaultSimulatorInterface::setSigma(sigma_ns);
103 
104 
106 
107  if (pmtParameters.getQE() != 1.0) {
108 
109  WARNING("Correct background rates with global efficiency " << pmtParameters.getQE() << endl);
110 
111  rates_Hz.correct(pmtParameters.getQE());
112  }
113 
114  DEBUG("PMT parameters: " << endl << pmtParameters << endl);
115  DEBUG("K40 rates: " << endl << rates_Hz << endl);
116 
118 
119  try {
120  load(detectorFile, detector);
121  }
122  catch(const JException& error) {
123  FATAL(error);
124  }
125 
126  JPMTParametersMap::Throw(false);
127 
129  JSummaryRouter summaryRouter;
130 
131  if (runbyrun.is_valid()) {
132 
133  NOTICE("Using run-by-run:" << endl << runbyrun << endl);
134 
135  if (!runbyrun.hasNext()) {
136  FATAL("Run-by-run simulation yields no input." << endl);
137  }
138 
139  if (rates_Hz.getSinglesRate() != 0.0) {
140  WARNING("Run-by-run simulation discards singles rate [Hz] " << rates_Hz.getSinglesRate() << endl);
141  }
142 
143  try {
144  simbad.reset(new JK40RunByRunSimulator(summaryRouter, rates_Hz));
145  simbad.reset(new JPMTDefaultSimulator(pmtParameters, detector));
146  simbad.reset(new JCLBDefaultSimulator());
147  }
148  catch(const JException& error) {
149  FATAL(error.what() << endl);
150  }
151 
152  } else {
153 
154  NOTICE("Using fixed rates [Hz]: " << rates_Hz << endl);
155 
156  try {
157  simbad.reset(new JK40DefaultSimulator(rates_Hz));
158  simbad.reset(fast ?
160  new JPMTDefaultSimulator(pmtParameters, detector));
161  simbad.reset(new JCLBDefaultSimulator());
162  }
163  catch(const JException& error) {
164  FATAL(error.what() << endl);
165  }
166  }
167 
168 
169  JTimer timerco("constructor");
170  JTimer timerrc("recycle");
171  JTimer timerIO("I/O");
172 
173  for (size_t i = 0; i != outputFile.size(); ++i) {
174 
175  outputFile[i].open();
176 
177  if (!outputFile[i].is_open()) {
178  FATAL("Error opening file " << outputFile[i] << endl);
179  }
180  }
181 
182  outputFile[0].put(JMeta(argc, argv));
183  outputFile[0].put(*gRandom);
184 
185  int counter = 0;
186 
187  for ( ; counter != numberOfSlices; ) {
188 
189  STATUS("slice: " << setw(10) << counter << '\r'); DEBUG(endl);
190 
191  int frame_index = counter + 1;
192 
193  if (runbyrun.hasNext()) {
194 
195  summaryRouter.update(runbyrun.next());
196 
197  summaryRouter.correct(dynamic_cast<const JPMTDefaultSimulatorInterface&>(simbad.getPMTSimulator()));
198 
199  frame_index = summaryRouter.getFrameIndex();
200  run = summaryRouter.getRunNumber();
201  }
202 
203  JDAQUTCExtended utc(JDAQUTCExtended::getInstance().getTimeNanoSecond() +
204  getTimeOfFrame(frame_index)); // ensure incremental UTC times (e.g. for JDAQSplit.cc)
205 
206  const JDAQChronometer chronometer(detector.getID(), run, frame_index, utc);
207 
208  timerco.start();
209 
210  JRandomTimeslice timeslice(chronometer, simbad);
211 
212  timerco.stop();
213 
214  timerIO.start();
215 
216  outputFile.rbegin()->put(timeslice); ++counter;
217 
218  timerIO.stop();
219 
220  for (size_t i = 1; i <= recycling.first && counter != numberOfSlices; ++i) {
221 
222  STATUS("slice: " << setw(10) << counter << '\r'); DEBUG(endl);
223 
224  timerrc.start();
225 
226  timeslice.recycle(recycling.second);
227 
228  timerrc.stop();
229 
230  timerIO.start();
231 
232  outputFile.rbegin()->put(timeslice); ++counter;
233 
234  timerIO.stop();
235  }
236  }
237  STATUS(endl);
238 
239  if (debug >= notice_t) {
240  timerco.print(cout, true);
241  timerrc.print(cout, true);
242  timerIO.print(cout, true);
243  }
244 
245  Head header;
246 
247  {
248  JHead buffer;
249 
250  buffer.DAQ.livetime_s = getLivetime(runbyrun->getFilelist());
251  buffer.K40.livetime_s = counter * getFrameTime() * 1.0e-9;
252 
253  buffer.push(&JHead::DAQ);
254  buffer.push(&JHead::K40);
255 
256  copy(buffer, header);
257  }
258 
259  outputFile[0].put(header);
260  outputFile[0].put(*gRandom);
261 
262  for (size_t i = 0; i != outputFile.size(); ++i) {
263  outputFile[i].close();
264  }
265 }
double getSigma(vector< double > &v)
get standard deviation of vector content
string outputFile
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62
#define STATUS(A)
Definition: JMessage.hh:63
#define NOTICE(A)
Definition: JMessage.hh:64
#define FATAL(A)
Definition: JMessage.hh:67
int debug
debug level
Definition: JSirene.cc:69
#define WARNING(A)
Definition: JMessage.hh:65
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:2158
Monte Carlo run header.
Definition: JHead.hh:1236
JAANET::DAQ DAQ
Definition: JHead.hh:1607
void push(T JHead::*pd)
Push given data member to Head.
Definition: JHead.hh:1374
JAANET::K40 K40
Definition: JHead.hh:1608
Detector data structure.
Definition: JDetector.hh:96
Default implementation of the simulation of K40 background.
Auxiliary class for map of PMT parameters.
double getQE(const JPMTIdentifier &id) const
Get QE of given PMT.
Auxiliary class for CPU timing and usage.
Definition: JTimer.hh:33
General exception.
Definition: JException.hh:24
virtual const char * what() const override
Get error message.
Definition: JException.hh:64
virtual const pointer_type & next() override
Get next element.
virtual bool hasNext() override
Check availability of next element.
Utility class to parse command line options.
Definition: JParser.hh:1714
Object writing to file.
Router for fast addressing of summary data in KM3NETDAQ::JDAQSummaryslice data structure as a functio...
void update(const JDAQSummaryslice *ps)
Update router.
K40 simulation based on run-by-run information.
unsigned int JTDC_t
leading edge [ns]
Definition: JDAQHit.hh:39
Data structure for UTC time.
void copy(const Head &from, JHead &to)
Copy header from from to to.
Definition: JHead.cc:162
void load(const std::string &file_name, JDetector &detector)
Load detector from input file.
@ debug_t
debug
Definition: JMessage.hh:29
@ notice_t
notice
Definition: JMessage.hh:32
T & getInstance(const T &object)
Get static instance from temporary object.
Definition: JObject.hh:75
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
double getLivetime(const std::string &file_name)
Get data taking live time.
KM3NeT DAQ data structures and auxiliaries.
Definition: DataQueue.cc:39
double getFrameTime()
Get frame time duration.
Definition: JDAQClock.hh:162
void setDAQLongprint(const bool option)
Set DAQ print option.
Definition: JDAQPrint.hh:28
double getTimeOfFrame(const int frame_index)
Get start time of frame in ns since start of run for a given frame index.
Definition: JDAQClock.hh:185
Definition: JSTDTypes.hh:14
The Head class reflects the header of Monte-Carlo event files, which consists of keys (also referred ...
Definition: Head.hh:65
double livetime_s
Live time [s].
Definition: JHead.hh:1053
double livetime_s
Live time [s].
Definition: JHead.hh:1107
Detector file.
Definition: JHead.hh:227
Empty structure for specification of parser element that is initialised (i.e. does not require input)...
Definition: JParser.hh:84
Auxiliary class for K40 rates.
Definition: JK40Rates.hh:41
double getSinglesRate() const
Get singles rate.
Definition: JK40Rates.hh:71
void correct(const double QE)
Correct rates for global efficiency,.
Definition: JK40Rates.hh:130
Auxiliary class for ROOT I/O of application specific meta data.
Definition: JMeta.hh:72
Auxiliary class to select summary data (KM3NETDAQ::JDAQSummaryslice) from the specified raw data file...
Definition: JRunByRun.hh:34
bool is_valid() const
Check validity of run by run options.
Definition: JRunByRun.hh:48
Timeslice with random data.