Jpp  15.0.0
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 "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 "km3net-dataformat/online/JDAQClock.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/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.

Author
mdejong

Definition in file software/JTimeslice/JRandomTimesliceWriter.cc.

Function Documentation

int main ( int  argc,
char **  argv 
)

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

40 {
41  using namespace std;
42  using namespace JPP;
43  using namespace KM3NETDAQ;
44 
46  string detectorFile;
47  Long64_t numberOfSlices;
48  JDAQHit::JTDC_t TCLB_ns;
49  JPMTParametersMap pmtParameters;
50  JK40Rates rates_Hz;
51  int run;
52  JRunByRun runbyrun;
53  pair<size_t, double> recycling;
54  UInt_t seed;
55  int debug;
56 
57  try {
58 
59  JParser<> zap("Auxiliary program to write time slices with random data.");
60 
61  zap['o'] = make_field(outputFile, "output file") = "timeslice.root";
62  zap['n'] = make_field(numberOfSlices);
63  zap['a'] = make_field(detectorFile, "detector.");
64  zap['R'] = make_field(run, "run number") = -1;
65  zap['r'] = make_field(runbyrun, "option for run-by-run mode") = JPARSER::initialised();
66  zap['P'] = make_field(pmtParameters, "PMT simulation data (or corresponding file name)") = JPARSER::initialised();
67  zap['B'] = make_field(rates_Hz, "background rates [Hz]") = JPARSER::initialised();
68  zap['T'] = make_field(TCLB_ns, "CLB state-machine time jitter") = 256; // [ns]
69  zap['N'] = make_field(recycling, "number of recycles / time interval for sampling data") = make_pair(0, 0.0);
70  zap['S'] = make_field(seed, "seed") = 0;
71  zap['d'] = make_field(debug, "debug") = 0;
72 
73  zap(argc, argv);
74  }
75  catch(const exception &error) {
76  FATAL(error.what() << endl);
77  }
78 
79 
80  gRandom->SetSeed(seed);
81 
82 
83  cout.tie(&cerr);
84 
86 
87  if (pmtParameters.getQE() != 1.0) {
88 
89  WARNING("Correct background rates with global efficiency " << pmtParameters.getQE() << endl);
90 
91  rates_Hz.correct(pmtParameters.getQE());
92  }
93 
94  DEBUG("PMT parameters: " << endl << pmtParameters << endl);
95  DEBUG("K40 rates: " << endl << rates_Hz << endl);
96 
98 
99  try {
100  load(detectorFile, detector);
101  }
102  catch(const JException& error) {
103  FATAL(error);
104  }
105 
106  JPMTParametersMap::Throw(false);
107 
109  JSummaryRouter summaryRouter;
110 
111  if (runbyrun.is_valid()) {
112 
113  NOTICE("Using run-by-run:" << endl << runbyrun << endl);
114 
115  if (!runbyrun.hasNext()) {
116  FATAL("Run-by-run simulation yields no input." << endl);
117  }
118 
119  if (rates_Hz.getSinglesRate() != 0.0) {
120  WARNING("Run-by-run simulation discards singles rate [Hz] " << rates_Hz.getSinglesRate() << endl);
121  }
122 
123  try {
124  simbad.reset(new JK40RunByRunSimulator(summaryRouter, rates_Hz));
125  simbad.reset(new JPMTRunByRunSimulator(summaryRouter, pmtParameters, detector));
126  simbad.reset(new JCLBDefaultSimulator());
127  }
128  catch(const JException& error) {
129  FATAL(error.what() << endl);
130  }
131 
132  } else {
133 
134  NOTICE("Using fixed rates [Hz]: " << rates_Hz << endl);
135 
136  try {
137  simbad.reset(new JK40DefaultSimulator(rates_Hz));
138  simbad.reset(new JPMTDefaultSimulator(pmtParameters, detector));
139  simbad.reset(new JCLBDefaultSimulator());
140  }
141  catch(const JException& error) {
142  FATAL(error.what() << endl);
143  }
144  }
145 
146 
147  JTimer timerco("constructor");
148  JTimer timerrc("recycle");
149  JTimer timerIO("I/O");
150 
151 
152  outputFile.open();
153 
154  if (!outputFile.is_open()) {
155  FATAL("Error opening file " << outputFile << endl);
156  }
157 
158  outputFile.put(*gRandom);
159  outputFile.put(JMeta(argc, argv));
160 
161  int counter = 0;
162 
163  for ( ; counter != numberOfSlices; ) {
164 
165  STATUS("slice: " << setw(10) << counter << '\r'); DEBUG(endl);
166 
167  int frame_index = counter + 1;
168 
169  if (runbyrun.hasNext()) {
170 
171  summaryRouter.update(runbyrun.next());
172 
173  summaryRouter.correct(dynamic_cast<const JPMTDefaultSimulatorInterface&>(simbad.getPMTSimulator()));
174 
175  frame_index = summaryRouter.getFrameIndex();
176  run = summaryRouter.getRunNumber();
177  }
178 
179  const JDAQChronometer chronometer(detector.getID(), run, frame_index, JDAQUTCExtended(getTimeOfFrame(frame_index)));
180 
181  timerco.start();
182 
183  JRandomTimeslice timeslice(chronometer, simbad);
184 
185  timerco.stop();
186 
187  timerIO.start();
188 
189  outputFile.put(timeslice); ++counter;
190 
191  timerIO.stop();
192 
193  for (size_t i = 1; i <= recycling.first && counter != numberOfSlices; ++i) {
194 
195  STATUS("slice: " << setw(10) << counter << '\r'); DEBUG(endl);
196 
197  timerrc.start();
198 
199  timeslice.recycle(recycling.second);
200 
201  timerrc.stop();
202 
203  outputFile.put(timeslice); ++counter;
204  }
205  }
206  STATUS(endl);
207 
208  if (debug >= notice_t) {
209  timerco.print(cout, true);
210  timerrc.print(cout, true);
211  timerIO.print(cout, true);
212  }
213 
214  outputFile.put(*gRandom);
215  outputFile.close();
216 }
Auxiliary class for ROOT I/O of application specific meta data.
Definition: JMeta.hh:70
Auxiliary class to select summary data (KM3NETDAQ::JDAQSummaryslice) from the specified raw data file...
Definition: JRunByRun.hh:32
Object writing to file.
Utility class to parse command line options.
Definition: JParser.hh:1500
General exception.
Definition: JException.hh:23
#define WARNING(A)
Definition: JMessage.hh:65
debug
Definition: JMessage.hh:29
Default implementation of the simulation of K40 background.
#define STATUS(A)
Definition: JMessage.hh:63
Detector data structure.
Definition: JDetector.hh:89
notice
Definition: JMessage.hh:32
Empty structure for specification of parser element that is initialised (i.e. does not require input)...
Definition: JParser.hh:66
string outputFile
Data structure for UTC time.
unsigned int JTDC_t
leading edge [ns]
Definition: JDAQHit.hh:39
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
K40 simulation based on run-by-run information.
void setDAQLongprint(const bool option)
Set DAQ print option.
Definition: JDAQPrint.hh:28
Detector file.
Definition: JHead.hh:196
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1961
Auxiliary class for CPU timing and usage.
Definition: JTimer.hh:32
#define NOTICE(A)
Definition: JMessage.hh:64
Auxiliary class for map of PMT parameters.
int debug
debug level
Definition: JSirene.cc:63
Router for fast addressing of summary data in JDAQSummaryslice data structure as a function of the op...
#define FATAL(A)
Definition: JMessage.hh:67
PMT simulation based on run-by-run information.
virtual const char * what() const override
Get error message.
Definition: JException.hh:48
void load(const std::string &file_name, JDetector &detector)
Load detector from input file.
then usage $script< input_file >< detector_file > fi set_variable OUTPUT_DIR set_variable SELECTOR JDAQTimesliceL1 set_variable DEBUG case set_variable DEBUG
do set_variable DETECTOR_TXT $WORKDIR detector
Timeslice with random data.
Auxiliary class for K40 rates.
Definition: JK40Rates.hh:41