Jpp  17.0.0-rc.1
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  double sigma_ns;
55  bool fast;
56  UInt_t seed;
57  int debug;
58 
59  try {
60 
61  JParser<> zap("Auxiliary program to write time slices with random data.");
62 
63  zap['o'] = make_field(outputFile, "output file") = "timeslice.root";
64  zap['n'] = make_field(numberOfSlices);
65  zap['a'] = make_field(detectorFile, "detector.");
66  zap['R'] = make_field(run, "run number") = -1;
67  zap['r'] = make_field(runbyrun, "option for run-by-run mode") = JPARSER::initialised();
68  zap['P'] = make_field(pmtParameters, "PMT simulation data (or corresponding file name)") = JPARSER::initialised();
69  zap['B'] = make_field(rates_Hz, "background rates [Hz]") = JPARSER::initialised();
70  zap['T'] = make_field(TCLB_ns, "CLB state-machine time jitter") = 256; // [ns]
71  zap['N'] = make_field(recycling, "number of recycles / time interval for sampling data") = make_pair(0, 0.0);
72  zap['s'] = make_field(sigma_ns, "intrinsic time smearing of K40 coincidences [ns]") = JK40DefaultSimulatorInterface::getSigma();
73  zap['F'] = make_field(fast, "fast - disable PMT simulation");
74  zap['S'] = make_field(seed, "seed") = 0;
75  zap['d'] = make_field(debug, "debug") = 0;
76 
77  zap(argc, argv);
78  }
79  catch(const exception &error) {
80  FATAL(error.what() << endl);
81  }
82 
83 
84  gRandom->SetSeed(seed);
85 
86  JK40DefaultSimulatorInterface::setSigma(sigma_ns);
87 
88  cout.tie(&cerr);
89 
91 
92  if (pmtParameters.getQE() != 1.0) {
93 
94  WARNING("Correct background rates with global efficiency " << pmtParameters.getQE() << endl);
95 
96  rates_Hz.correct(pmtParameters.getQE());
97  }
98 
99  DEBUG("PMT parameters: " << endl << pmtParameters << endl);
100  DEBUG("K40 rates: " << endl << rates_Hz << endl);
101 
103 
104  try {
105  load(detectorFile, detector);
106  }
107  catch(const JException& error) {
108  FATAL(error);
109  }
110 
111  JPMTParametersMap::Throw(false);
112 
114  JSummaryRouter summaryRouter;
115 
116  if (runbyrun.is_valid()) {
117 
118  NOTICE("Using run-by-run:" << endl << runbyrun << endl);
119 
120  if (!runbyrun.hasNext()) {
121  FATAL("Run-by-run simulation yields no input." << endl);
122  }
123 
124  if (rates_Hz.getSinglesRate() != 0.0) {
125  WARNING("Run-by-run simulation discards singles rate [Hz] " << rates_Hz.getSinglesRate() << endl);
126  }
127 
128  try {
129  simbad.reset(new JK40RunByRunSimulator(summaryRouter, rates_Hz));
130  simbad.reset(new JPMTRunByRunSimulator(summaryRouter, pmtParameters, detector));
131  simbad.reset(new JCLBDefaultSimulator());
132  }
133  catch(const JException& error) {
134  FATAL(error.what() << endl);
135  }
136 
137  } else {
138 
139  NOTICE("Using fixed rates [Hz]: " << rates_Hz << endl);
140 
141  try {
142  simbad.reset(new JK40DefaultSimulator(rates_Hz));
143  simbad.reset(fast ?
145  new JPMTDefaultSimulator(pmtParameters, detector));
146  simbad.reset(new JCLBDefaultSimulator());
147  }
148  catch(const JException& error) {
149  FATAL(error.what() << endl);
150  }
151  }
152 
153 
154  JTimer timerco("constructor");
155  JTimer timerrc("recycle");
156  JTimer timerIO("I/O");
157 
158 
159  outputFile.open();
160 
161  if (!outputFile.is_open()) {
162  FATAL("Error opening file " << outputFile << endl);
163  }
164 
165  outputFile.put(*gRandom);
166  outputFile.put(JMeta(argc, argv));
167 
168  int counter = 0;
169 
170  for ( ; counter != numberOfSlices; ) {
171 
172  STATUS("slice: " << setw(10) << counter << '\r'); DEBUG(endl);
173 
174  int frame_index = counter + 1;
175 
176  if (runbyrun.hasNext()) {
177 
178  summaryRouter.update(runbyrun.next());
179 
180  summaryRouter.correct(dynamic_cast<const JPMTDefaultSimulatorInterface&>(simbad.getPMTSimulator()));
181 
182  frame_index = summaryRouter.getFrameIndex();
183  run = summaryRouter.getRunNumber();
184  }
185 
186  const JDAQChronometer chronometer(detector.getID(), run, frame_index, JDAQUTCExtended(getTimeOfFrame(frame_index)));
187 
188  timerco.start();
189 
190  JRandomTimeslice timeslice(chronometer, simbad);
191 
192  timerco.stop();
193 
194  timerIO.start();
195 
196  outputFile.put(timeslice); ++counter;
197 
198  timerIO.stop();
199 
200  for (size_t i = 1; i <= recycling.first && counter != numberOfSlices; ++i) {
201 
202  STATUS("slice: " << setw(10) << counter << '\r'); DEBUG(endl);
203 
204  timerrc.start();
205 
206  timeslice.recycle(recycling.second);
207 
208  timerrc.stop();
209 
210  timerIO.start();
211 
212  outputFile.put(timeslice); ++counter;
213 
214  timerIO.stop();
215  }
216  }
217  STATUS(endl);
218 
219  if (debug >= notice_t) {
220  timerco.print(cout, true);
221  timerrc.print(cout, true);
222  timerIO.print(cout, true);
223  }
224 
225  outputFile.put(*gRandom);
226  outputFile.close();
227 }
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
double getSigma(vector< double > &v)
get standard deviation of vector content
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:224
#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:66
Router for fast addressing of summary data in KM3NETDAQ::JDAQSummaryslice data structure as a functio...
#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.
do set_variable DETECTOR_TXT $WORKDIR detector
Timeslice with random data.
Auxiliary class for K40 rates.
Definition: JK40Rates.hh:41
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62