Jpp  17.2.1-pre0
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
JTriggerEfficiency.cc
Go to the documentation of this file.
1 
2 #include <string>
3 #include <iostream>
4 #include <iomanip>
5 #include <vector>
6 #include <set>
7 
8 #include "TH1D.h"
9 #include "TRandom3.h"
10 
15 
16 #include "JDAQ/JDAQTimesliceIO.hh"
17 #include "JDAQ/JDAQEventIO.hh"
19 #include "JDAQ/JDAQToolkit.hh"
22 
23 #include "JAAnet/JHead.hh"
24 #include "JAAnet/JHeadToolkit.hh"
25 #include "JAAnet/JAAnetToolkit.hh"
26 
27 #include "JPhysics/JK40Rates.hh"
28 
29 #include "JDetector/JDetector.hh"
33 #include "JDetector/JPMTRouter.hh"
34 #include "JDetector/JTimeRange.hh"
39 
40 #include "JTrigger/JHit.hh"
41 #include "JTrigger/JHitToolkit.hh"
42 #include "JTrigger/JTimeslice.hh"
45 #include "JTrigger/JHitL0.hh"
46 #include "JTrigger/JHitL1.hh"
47 #include "JTrigger/JBuildL1.hh"
48 #include "JTrigger/JBuildL2.hh"
52 #include "JTrigger/JTriggerBits.hh"
56 #include "JTrigger/JTimesliceL1.hh"
64 
70 #include "JSupport/JSupport.hh"
71 #include "JSupport/JRunByRun.hh"
72 #include "JSupport/JMeta.hh"
73 
74 #include "Jeep/JPrint.hh"
75 #include "Jeep/JParser.hh"
76 #include "Jeep/JMessage.hh"
77 
78 
79 /**
80  * \file
81  * Auxiliary program to trigger Monte Carlo events.
82  *
83  * The options
84  * - <tt>-a</tt> refers to the detector configuration that is used to convert Monte Carlo true information to raw data; and
85  * - <tt>-b</tt> refers to the detector configuration that is used to convert raw data to calibrated data.
86  *
87  * The event counter of the triggered events, as returned by KM3NETDAQ::JDAQEvent::getCounter(),
88  * is set to the corresponding index of the Monte Carlo event in the output.\n
89  * This allows for correlating the DAQ event to the Monte Carlo event.\n
90  * The time of the DAQ hits can be correlated to the time of the Monte Carlo hits using the frame index
91  * and the time of the Monte Carlo event (Evt::mc_t), respectively (see also time_converter).\n
92  * The universal time of the DAQ event is set to the universal time of the Monte Carlo event
93  * with a correction for the time of the event (read hits) within the time slice.\n
94  * In run-by-run mode, the trigger parameters as well as the singles rates and status of the PMTs
95  * are read from the given raw data file.\n
96  * The two-fold (and higher) coincidence rates should still be provided on the command line.
97  *
98  * \author mdejong
99  */
100 int main(int argc, char **argv)
101 {
102  using namespace std;
103  using namespace JPP;
104  using namespace KM3NETDAQ;
105 
107 
108  JMultipleFileScanner<> inputFile;
110  JLimit_t& numberOfEvents = inputFile.getLimit();
111  string detectorFileA;
112  string detectorFileB;
113  int run;
115  bool triggeredEventsOnly;
116  JPMTParametersMap pmtParameters;
117  JK40Rates rates_Hz;
118  JRunByRun runbyrun;
119  double sigma_ns;
120  UInt_t seed;
121  int debug;
122 
123  try {
124 
125  JParser<> zap("Auxiliary program to trigger Monte Carlo events.");
126 
127  zap['f'] = make_field(inputFile, "input file (output of detector simulation)");
128  zap['o'] = make_field(outputFile, "output file") = "trigger_efficieny.root";
129  zap['n'] = make_field(numberOfEvents) = JLimit::max();
130  zap['a'] = make_field(detectorFileA, "detector used for conversion from Monte Carlo truth to raw data.");
131  zap['b'] = make_field(detectorFileB, "detector used for conversion of raw data to calibrated data.") = "";
132  zap['R'] = make_field(run, "run number") = -1;
133  zap['r'] = make_field(runbyrun, "option for run-by-run mode") = JPARSER::initialised();
134  zap['@'] = make_field(parameters, "Trigger parameters (or corresponding file name)") = JPARSER::initialised();
135  zap['O'] = make_field(triggeredEventsOnly, "optionally write only triggered events.");
136  zap['P'] = make_field(pmtParameters, "PMT simulation data (or corresponding file name)") = JPARSER::initialised();
137  zap['B'] = make_field(rates_Hz, "background rates [Hz]") = JPARSER::initialised();
138  zap['s'] = make_field(sigma_ns, "intrinsic time smearing of K40 coincidences [ns]") = JK40DefaultSimulatorInterface::getSigma();
139  zap['S'] = make_field(seed, "seed") = 0;
140  zap['d'] = make_field(debug, "debug") = 0;
141 
142  zap(argc, argv);
143  }
144  catch(const exception &error) {
145  FATAL(error.what() << endl);
146  }
147 
148  gRandom->SetSeed(seed);
149 
150  JK40DefaultSimulatorInterface::setSigma(sigma_ns);
151 
152  cout.tie(&cerr);
153 
155 
156  if (detectorFileB == "") {
157  detectorFileB = detectorFileA;
158  }
159 
160 
161  JDetector detectorA;
162  JDetector detectorB;
163 
164  try {
165  load(detectorFileA, detectorA);
166  load(detectorFileB, detectorB);
167  }
168  catch(const JException& error) {
169  FATAL(error);
170  }
171 
172  JPMTParametersMap::Throw(true);
173 
174  if (!pmtParameters.is_valid()) {
175  FATAL("Invalid PMT parameters " << pmtParameters << endl);
176  }
177 
178  if (pmtParameters.getQE() != 1.0) {
179 
180  WARNING("Correct background rates with global efficiency " << pmtParameters.getQE() << endl);
181 
182  rates_Hz.correct(pmtParameters.getQE());
183  }
184 
185  Head header;
186 
187  try {
188  header = getHeader(inputFile);
189  }
190  catch(const JException& error) {
191  FATAL(error);
192  }
193 
194  const JModuleRouter moduleRouter(detectorB);
195  JDetectorSimulator simbad (detectorA);
196  JSummaryRouter summaryRouter;
197 
198  if (runbyrun.is_valid()) {
199 
200  NOTICE("Using run-by-run:" << endl << runbyrun << endl);
201 
202  if (!runbyrun.hasNext()) {
203  FATAL("Run-by-run simulation yields no input." << endl);
204  }
205 
206  if (rates_Hz.getSinglesRate() != 0.0) {
207  WARNING("Run-by-run simulation discards singles rate [Hz] " << rates_Hz.getSinglesRate() << endl);
208  }
209 
210  try {
211  simbad.reset(new JK40RunByRunSimulator(summaryRouter, rates_Hz));
212  simbad.reset(new JPMTRunByRunSimulator(summaryRouter, pmtParameters, detectorA));
213  simbad.reset(new JCLBRunByRunSimulator(summaryRouter));
214  }
215  catch(const JException& error) {
216  FATAL(error.what() << endl);
217  }
218 
219  try {
220 
221  parameters = getTriggerParameters(runbyrun->getFilelist());
222 
223  NOTICE("Set trigger parameters from run-by-run input." << endl);
224  }
225  catch(const JException& error) {
226  WARNING("No trigger parameters from run-by-run input;\nrun with default/user input." << endl);
227  }
228 
229  // set live time
230 
231  JHead buffer(header);
232 
233  buffer.DAQ.livetime_s = getLivetime(runbyrun->getFilelist());
234  buffer.push(&JHead::DAQ);
235 
236  copy(buffer, header);
237 
238  } else {
239 
240  NOTICE("Using fixed rates [Hz]: " << rates_Hz << endl);
241 
242  try {
243  simbad.reset(new JK40DefaultSimulator(rates_Hz));
244  simbad.reset(new JPMTDefaultSimulator(pmtParameters, detectorA));
245  simbad.reset(new JCLBDefaultSimulator());
246  }
247  catch(const JException& error) {
248  FATAL(error.what() << endl);
249  }
250  }
251 
252  // detector
253 
254  if (parameters.disableHighRateVeto) {
255 
256  NOTICE("Disabling high-rate veto of all PMTs." << endl);
257 
258  detectorB.setPMTStatus(HIGH_RATE_VETO_DISABLE);
259  }
260 
261  parameters.set(getMaximalDistance(detectorB));
262 
263  DEBUG("Trigger:" << endl << parameters << endl);
264  DEBUG("PMT parameters:" << endl << pmtParameters << endl);
265 
266  const double Tmax = max(getMaximalTime(detectorA),
267  getMaximalTime(detectorB));
268 
269  const JTimeRange period(-(Tmax + parameters.TMaxLocal_ns),
270  +(Tmax + parameters.TMaxLocal_ns));
271 
272  typedef double hit_type;
273 
274  typedef JSuperFrame1D<hit_type> JSuperFrame1D_t;
275  typedef JSuperFrame2D<hit_type> JSuperFrame2D_t;
276  typedef JTimeslice <hit_type> JTimeslice_t;
277  typedef JBuildL1 <hit_type> JBuildL1_t;
278  typedef JBuildL2 <hit_type> JBuildL2_t;
279 
280  const JBuildL1_t buildL1(parameters);
281  const JBuildL2_t buildL2(parameters.L2);
282  const JBuildL2_t buildSN(parameters.SN);
283 
284  JTimesliceRouter timesliceRouter(parameters.numberOfBins);
285 
286  const JTrigger3DMuon trigger3DMuon (parameters);
287  const JTrigger3DShower trigger3DShower(parameters);
288  const JTriggerMXShower triggerMXShower(parameters, detectorB);
289 
290  const JPosition3D center = get<JPosition3D>(header);
291 
292  NOTICE("Apply detector offset from Monte Carlo run header (" << center << ")" << endl);
293 
294  detectorA -= center;
295  detectorB -= center;
296 
297  TH1D h1("Trigger bits", NULL, NUMBER_OF_TRIGGER_BITS, -0.5, NUMBER_OF_TRIGGER_BITS - 0.5);
298 
299 
300  outputFile.open();
301 
302  if (!outputFile.is_open()) {
303  FATAL("Error opening file " << outputFile << endl);
304  }
305 
306  outputFile.put(*gRandom);
307  outputFile.put(JMeta(argc, argv));
308  outputFile.put(header);
309  outputFile.put(parameters);
310 
311  JLimit_t limit = inputFile.getLimit();
312  counter_type number_of_events = 0;
313  int trigger_counter = 0;
314 
315  for (JMultipleFileScanner<>::const_iterator file = inputFile.begin(); file != inputFile.end(); ++file) {
316 
317  int mc_run_id = 0;
318 
319  try {
320 
321  const JHead head(getHeader(*file));
322 
323  mc_run_id = head.start_run.run_id;
324  }
325  catch(const JException& error) {
326  FATAL(error);
327  }
328 
330 
331  limit.setLowerLimit(limit.getLowerLimit() - in.skip(limit.getLowerLimit()));
332 
333  for ( ; in.hasNext() && number_of_events != limit; ++number_of_events) {
334 
335  STATUS("event: " << setw(10) << number_of_events << '\r'); DEBUG(endl);
336 
337  Evt* event = in.next();
338 
339  event->mc_run_id = mc_run_id;
340 
341  DEBUG(*event << endl);
342 
343  bool trigger = false;
344 
345  if (!event->mc_hits.empty()) {
346 
347  int frame_index = (int) in.getCounter() + 1; // 1 event / slice
348 
349  if (runbyrun.is_valid() && runbyrun.hasNext()) {
350 
351  summaryRouter.update(runbyrun.next());
352 
353  summaryRouter.correct(dynamic_cast<const JPMTDefaultSimulatorInterface&>(simbad.getPMTSimulator()));
354 
355  frame_index = summaryRouter.getFrameIndex();
356  run = summaryRouter.getRunNumber();
357  }
358 
359  // set the event time!
360 
361  JTimeRange timeRange = getTimeRange(*event, period);
362 
363  if (!timeRange.is_valid()) {
364  timeRange.setRange(0.0,0.0);
365  }
366 
367  const double t0 = 0.5 * (timeRange.getLowerLimit() + timeRange.getUpperLimit());
368  const double t1 = gRandom->Rndm() * getFrameTime();
369 
370  event->mc_t = getTimeOfFrame(frame_index) + t1 - t0; // time since start of data taking run
371 
372  timeRange.add(event->mc_t);
373  timeRange.add(period);
374 
375  JDAQUTCExtended utc(getTimeOfFrame(trigger_counter + 1)); // ensure incremental UTC times (e.g. for JDAQSplit.cc)
376 
377  if (event->mc_event_time != TTimeStamp(0)) {
378  utc = getDAQUTCExtended(event->mc_event_time, t1); // set UTC time according Monte Carlo event
379  }
380 
381  const JDAQChronometer chronometer(detectorB.getID(), (run != -1 ? run : mc_run_id), frame_index, utc);
382 
383  const JEventTimeslice timeslice(chronometer, simbad, *event, period);
384 
385  DEBUG(timeslice << endl);
386 
387 
388  timesliceRouter.configure(timeslice);
389 
390  JTimeslice_t timesliceL0(timeslice.getDAQChronometer());
391  JTimeslice_t timesliceL1(timeslice.getDAQChronometer());
392  JTimeslice_t timesliceL2(timeslice.getDAQChronometer());
393  JTimeslice_t timesliceSN(timeslice.getDAQChronometer());
394 
395  for (JDAQTimeslice::const_iterator super_frame = timeslice.begin(); super_frame != timeslice.end(); ++super_frame) {
396 
397  if (moduleRouter.hasModule(super_frame->getModuleID())) {
398 
399  // calibration
400 
401  const JModule& module = moduleRouter.getModule(super_frame->getModuleID());
402  const JSuperFrame2D_t& buffer = JSuperFrame2D_t::demultiplex(*super_frame, module);
403 
404  // L0
405 
406  timesliceL0.push_back(JSuperFrame1D_t(buffer));
407 
408  // L1
409 
410  timesliceL1.push_back(JSuperFrame1D_t(super_frame->getDAQChronometer(),
411  super_frame->getModuleIdentifier(),
412  module.getPosition()));
413 
414  buildL1(*timesliceL0.rbegin(), back_inserter(*timesliceL1.rbegin()));
415 
416  // L2
417 
418  timesliceL2.push_back(JSuperFrame1D_t(super_frame->getDAQChronometer(),
419  super_frame->getModuleIdentifier(),
420  module.getPosition()));
421 
422  buildL2(buffer, *timesliceL1.rbegin(), back_inserter(*timesliceL2.rbegin()));
423 
424  // SN
425 
426  timesliceSN.push_back(JSuperFrame1D_t(super_frame->getDAQChronometer(),
427  super_frame->getModuleIdentifier(),
428  module.getPosition()));
429 
430  buildSN(buffer, *timesliceL1.rbegin(), back_inserter(*timesliceSN.rbegin()));
431 
432  DEBUG("L0 " << setw(8) << timesliceL0.rbegin()->getModuleID() << ' ' << setw(8) << timesliceL0.rbegin()->size() << endl);
433  DEBUG("L1 " << setw(8) << timesliceL1.rbegin()->getModuleID() << ' ' << setw(8) << timesliceL1.rbegin()->size() << endl);
434  DEBUG("L2 " << setw(8) << timesliceL2.rbegin()->getModuleID() << ' ' << setw(8) << timesliceL2.rbegin()->size() << endl);
435  DEBUG("SN " << setw(8) << timesliceSN.rbegin()->getModuleID() << ' ' << setw(8) << timesliceSN.rbegin()->size() << endl);
436  }
437  }
438 
439 
440  // Trigger
441 
442  JTriggerInput trigger_input(timesliceL2);
443  JTriggerOutput trigger_output;
444 
445  trigger3DMuon (trigger_input, back_inserter(trigger_output));
446  trigger3DShower(trigger_input, back_inserter(trigger_output));
447  triggerMXShower(trigger_input, timesliceL0, back_inserter(trigger_output));
448 
449  trigger_output.merge(JEventOverlap(parameters.TMaxEvent_ns));
450 
451  for (JTriggerOutput::const_iterator to = trigger_output.begin(); to != trigger_output.end(); ++to) {
452 
453  for (int i = 0; i != h1.GetNbinsX(); ++i) {
454  if (to->hasTriggerBit(i)) {
455  h1.Fill((double) i);
456  }
457  }
458 
459  JTimeRange eventTime = getTimeRange(*to).add(getTimeOfRTS(*to));
460 
461  DEBUG("Event time: "
462  << to->getFrameIndex() << ' '
463  << eventTime << ' '
464  << timeRange << endl);
465 
466  if (timeRange.overlap(eventTime)) {
467 
468  JTriggeredEvent tev(*to,
469  timesliceRouter,
470  moduleRouter,
471  parameters.TMaxLocal_ns,
473 
474  // Set the event counter to the index of the Monte Carlo event in the output.
475 
476  tev.setCounter(trigger_counter);
477 
478  outputFile.put(tev);
479 
480  trigger = true;
481  }
482  }
483 
484 
485  if (!triggeredEventsOnly || trigger) {
486 
487  if (parameters.writeL0()) {
488  outputFile.put(timeslice);
489  }
490 
491  if (parameters.writeL1()) {
492  outputFile.put(JTimesliceL1<JDAQTimesliceL1>(timesliceL1, timesliceRouter, moduleRouter, parameters.TMaxLocal_ns));
493  }
494 
495  if (parameters.writeL2()) {
496  outputFile.put(JTimesliceL1<JDAQTimesliceL2>(timesliceL2, timesliceRouter, moduleRouter, parameters.L2.TMaxLocal_ns));
497  }
498 
499  if (parameters.writeSN()) {
500  outputFile.put(JTimesliceL1<JDAQTimesliceSN>(timesliceSN, timesliceRouter, moduleRouter, parameters.SN.TMaxLocal_ns));
501  }
502 
503  if (parameters.writeSummary()) {
504  outputFile.put(JSummaryslice(chronometer,simbad));
505  }
506  }
507  }
508 
509  if (!triggeredEventsOnly || trigger) {
510 
511  outputFile.put(*event);
512 
513  ++trigger_counter;
514  }
515  }
516  }
517  STATUS(endl);
518 
519 
521 
522  io >> outputFile;
523 
524  outputFile.put(h1);
525  outputFile.put(*gRandom);
526  outputFile.close();
527 }
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
#define WARNING(A)
Definition: JMessage.hh:65
debug
Definition: JMessage.hh:29
T getLowerLimit() const
Get lower limit.
Definition: JRange.hh:202
int main(int argc, char *argv[])
Definition: Main.cc:15
ROOT TTree parameter settings of various packages.
static const unsigned int NUMBER_OF_TRIGGER_BITS
Number of trigger bits.
Default implementation of the simulation of K40 background.
const JModule & getModule(const JObjectID &id) const
Get module parameters.
Data structure for a composite optical module.
Definition: JModule.hh:68
Router for fast addressing of hits in KM3NETDAQ::JDAQTimeslice data structure as a function of the op...
#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
void configure(const JDAQTimeslice &timeslice)
Configure.
Recording of objects on file according a format that follows from the file name extension.
Router for direct addressing of module data in detector data structure.
static const int HIGH_RATE_VETO_DISABLE
Enable (disable) use of high-rate veto test if this status bit is 0 (1);.
Definition: pmt_status.hh:13
*fatal Wrong number of arguments esac JCookie sh typeset Z DETECTOR typeset Z SOURCE_RUN typeset Z TARGET_RUN set_variable PARAMETERS_FILE $WORKDIR parameters
Definition: diff-Tuna.sh:38
Empty structure for specification of parser element that is initialised (i.e. does not require input)...
Definition: JParser.hh:83
Long64_t counter_type
Type definition for counter.
double livetime_s
Live time [s].
Definition: JHead.hh:1040
const JPMTSimulator & getPMTSimulator() const
Get PMT simulator.
double getTimeOfRTS(const JDAQChronometer &chronometer)
Get time of last RTS in ns since start of run for a given chronometer.
Basic data structure for time and time over threshold information of hit.
string outputFile
Data structure for UTC time.
Data structure for detector geometry and calibration.
Tools for handling different hit types.
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
1-dimensional frame with time calibrated data from one optical module.
JTimeRange getTimeRange(const Evt &event)
Get time range (i.e. time between earliest and latest hit) of Monte Carlo event.
K40 simulation based on run-by-run information.
Head getHeader(const JMultipleFileScanner_t &file_list)
Get Monte Carlo header.
Basic data structure for L0 hit.
Type list.
Definition: JTypeList.hh:22
void setDAQLongprint(const bool option)
Set DAQ print option.
Definition: JDAQPrint.hh:28
Auxiliary class for defining the range of iterations of objects.
Definition: JLimit.hh:41
Template L2 builder.
Definition: JBuildL2.hh:45
void push(T JHead::*pd)
Push given data member to Head.
Definition: JHead.hh:1306
double getMaximalDistance(const JDetector &detector)
Get maximal distance between modules in detector.
I/O formatting auxiliaries.
Auxiliaries for creation of summary data.
virtual bool hasNext() override
Check availability of next element.
void merge(const JMatch_t &match)
Merge events.
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1993
CLB simulation based on run-by-run information.
range_type & add(argument_type x)
Add offset.
Definition: JRange.hh:446
double getFrameTime()
Get frame time duration.
Definition: JDAQClock.hh:162
void setRange(const range_type &range)
Set range.
Definition: JRange.hh:146
JAANET::start_run start_run
Definition: JHead.hh:1509
ROOT I/O of application specific meta data.
T getUpperLimit() const
Get upper limit.
Definition: JRange.hh:213
Definition of hit and track types and auxiliary methods for handling Monte Carlo data.
#define NOTICE(A)
Definition: JMessage.hh:64
Auxiliary class to create summary data.
const JDAQChronometer & getDAQChronometer() const
Get DAQ chronometer.
Support methods.
Auxiliary class to build JDAQEvent for a triggered event.
Q DAQ
Definition: JDataQuality.sh:53
Auxiliary class for map of PMT parameters.
Direct access to PMT in detector data structure.
const JPosition3D & getPosition() const
Get position.
Definition: JPosition3D.hh:130
double getMaximalTime(const double R_Hz)
Get maximal time for given rate.
Router for fast addressing of summary data in KM3NETDAQ::JDAQSummaryslice data structure as a functio...
General purpose messaging.
Monte Carlo run header.
Definition: JHead.hh:1167
int run_id
MC run number.
Definition: JHead.hh:143
The Head class reflects the header of Monte-Carlo event files, which consists of keys (also referred ...
Definition: Head.hh:65
Template L1 hit builder.
Definition: JBuildL1.hh:85
#define FATAL(A)
Definition: JMessage.hh:67
Scanning of objects from multiple files according a format that follows from the extension of each fi...
void reset(JK40Simulator *k40Simulator)
Reset K40 simulator.
PMT simulation based on run-by-run information.
void setCounter(const JTriggerCounter_t counter)
Set trigger counter.
bool is_valid() const
Check validity of range.
Definition: JRange.hh:311
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.
Auxiliary class to build JDAQTimeslice for L1 timeslice.
Definition: JTimesliceL1.hh:36
double getLivetime(const std::string &file_name)
Get data taking live time.
Utility class to parse command line options.
Map of associated modules in detector.
virtual const pointer_type & next() override
Get next element.
bool hasModule(const JObjectID &id) const
Has module.
void copy(const Head &from, JHead &to)
Copy header from from to to.
Definition: JHead.cc:161
2-dimensional frame with time calibrated data from one optical module.
Data structure for input to trigger algorithm.
counter_type getCounter() const
Get counter.
Data structure for position in three dimensions.
Definition: JPosition3D.hh:36
const JLimit & getLimit() const
Get limit.
Definition: JLimit.hh:73
Setting of trigger bits.
JAANET::DAQ DAQ
Definition: JHead.hh:1532
then fatal Wrong number of arguments fi set_variable DETECTOR $argv[1] set_variable INPUT_FILE $argv[2] eval JPrintDetector a $DETECTOR O IDENTIFIER eval JPrintDetector a $DETECTOR O SUMMARY JAcoustics sh $DETECTOR_ID source JAcousticsToolkit sh CHECK_EXIT_CODE typeset A EMITTERS get_tripods $WORKDIR tripod txt EMITTERS get_transmitters $WORKDIR transmitter txt EMITTERS for EMITTER in
Definition: JCanberra.sh:46
JDAQUTCExtended getDAQUTCExtended(const TTimeStamp &t0, const double t1=0.0)
Get DAQ UTC time.
Definition: JDAQToolkit.hh:26
Match of two events considering overlap in time.
bool overlap(const range_type &range) const
Test overlap with given range.
Definition: JRange.hh:382
Timeslice with Monte Carlo event.
JTriggerParameters getTriggerParameters(const JMultipleFileScanner_t &file_list)
Get trigger parameters.
Basic data structure for L1 hit.
int debug
debug level
Auxiliary class for K40 rates.
Definition: JK40Rates.hh:41
The Evt class respresent a Monte Carlo (MC) event as well as an offline event.
Definition: Evt.hh:20
virtual skip_type skip(const skip_type ns) override
Skip items.
Time slice with calibrated data.
Definition: JTimeslice.hh:26
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62