Jpp  17.3.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 [ns]") = 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 
90 
91  if (pmtParameters.getQE() != 1.0) {
92 
93  WARNING("Correct background rates with global efficiency " << pmtParameters.getQE() << endl);
94 
95  rates_Hz.correct(pmtParameters.getQE());
96  }
97 
98  DEBUG("PMT parameters: " << endl << pmtParameters << endl);
99  DEBUG("K40 rates: " << endl << rates_Hz << endl);
100 
102 
103  try {
104  load(detectorFile, detector);
105  }
106  catch(const JException& error) {
107  FATAL(error);
108  }
109 
110  JPMTParametersMap::Throw(false);
111 
113  JSummaryRouter summaryRouter;
114 
115  if (runbyrun.is_valid()) {
116 
117  NOTICE("Using run-by-run:" << endl << runbyrun << endl);
118 
119  if (!runbyrun.hasNext()) {
120  FATAL("Run-by-run simulation yields no input." << endl);
121  }
122 
123  if (rates_Hz.getSinglesRate() != 0.0) {
124  WARNING("Run-by-run simulation discards singles rate [Hz] " << rates_Hz.getSinglesRate() << endl);
125  }
126 
127  try {
128  simbad.reset(new JK40RunByRunSimulator(summaryRouter, rates_Hz));
129  simbad.reset(new JPMTRunByRunSimulator(summaryRouter, pmtParameters, detector));
130  simbad.reset(new JCLBDefaultSimulator());
131  }
132  catch(const JException& error) {
133  FATAL(error.what() << endl);
134  }
135 
136  } else {
137 
138  NOTICE("Using fixed rates [Hz]: " << rates_Hz << endl);
139 
140  try {
141  simbad.reset(new JK40DefaultSimulator(rates_Hz));
142  simbad.reset(fast ?
144  new JPMTDefaultSimulator(pmtParameters, detector));
145  simbad.reset(new JCLBDefaultSimulator());
146  }
147  catch(const JException& error) {
148  FATAL(error.what() << endl);
149  }
150  }
151 
152 
153  JTimer timerco("constructor");
154  JTimer timerrc("recycle");
155  JTimer timerIO("I/O");
156 
157 
158  outputFile.open();
159 
160  if (!outputFile.is_open()) {
161  FATAL("Error opening file " << outputFile << endl);
162  }
163 
164  outputFile.put(*gRandom);
165  outputFile.put(JMeta(argc, argv));
166 
167  int counter = 0;
168 
169  for ( ; counter != numberOfSlices; ) {
170 
171  STATUS("slice: " << setw(10) << counter << '\r'); DEBUG(endl);
172 
173  int frame_index = counter + 1;
174 
175  if (runbyrun.hasNext()) {
176 
177  summaryRouter.update(runbyrun.next());
178 
179  summaryRouter.correct(dynamic_cast<const JPMTDefaultSimulatorInterface&>(simbad.getPMTSimulator()));
180 
181  frame_index = summaryRouter.getFrameIndex();
182  run = summaryRouter.getRunNumber();
183  }
184 
185  const JDAQChronometer chronometer(detector.getID(), run, frame_index, JDAQUTCExtended(getTimeOfFrame(frame_index)));
186 
187  timerco.start();
188 
189  JRandomTimeslice timeslice(chronometer, simbad);
190 
191  timerco.stop();
192 
193  timerIO.start();
194 
195  outputFile.put(timeslice); ++counter;
196 
197  timerIO.stop();
198 
199  for (size_t i = 1; i <= recycling.first && counter != numberOfSlices; ++i) {
200 
201  STATUS("slice: " << setw(10) << counter << '\r'); DEBUG(endl);
202 
203  timerrc.start();
204 
205  timeslice.recycle(recycling.second);
206 
207  timerrc.stop();
208 
209  timerIO.start();
210 
211  outputFile.put(timeslice); ++counter;
212 
213  timerIO.stop();
214  }
215  }
216  STATUS(endl);
217 
218  if (debug >= notice_t) {
219  timerco.print(cout, true);
220  timerrc.print(cout, true);
221  timerIO.print(cout, true);
222  }
223 
224  outputFile.put(*gRandom);
225  outputFile.close();
226 }
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:1517
General exception.
Definition: JException.hh:23
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:83
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:226
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1993
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.
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
then echo WARNING
Definition: JTuneHV.sh:91
Timeslice with random data.
int debug
debug level
Auxiliary class for K40 rates.
Definition: JK40Rates.hh:41
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62