Jpp
Functions
JTimeConverter.cc File Reference
#include <string>
#include <iostream>
#include <iomanip>
#include <map>
#include <set>
#include "km3net-dataformat/offline/Head.hh"
#include "km3net-dataformat/offline/Evt.hh"
#include "km3net-dataformat/offline/Hit.hh"
#include "JDAQ/JDAQEventIO.hh"
#include "JAAnet/JAAnetToolkit.hh"
#include "JDetector/JDetector.hh"
#include "JDetector/JDetectorToolkit.hh"
#include "JDetector/JModuleRouter.hh"
#include "JTrigger/JTimeConverter.hh"
#include "JSupport/JTriggeredFileScanner.hh"
#include "JSupport/JMonteCarloFileSupportkit.hh"
#include "JSupport/JSupport.hh"
#include "Jeep/JPrint.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

Example program to test conversion between Monte Carlo and DAQ times.

Author
mdejong

Definition in file JTimeConverter.cc.

Function Documentation

◆ main()

int main ( int  argc,
char **  argv 
)

Definition at line 70 of file JTimeConverter.cc.

71 {
72  using namespace std;
73  using namespace JPP;
74  using namespace KM3NETDAQ;
75 
76  JTriggeredFileScanner<> inputFile;
77  JLimit_t& numberOfEvents = inputFile.getLimit();
78  string detectorFile;
79  double Tmax_ns;
80  int debug;
81 
82  try {
83 
84  JParser<> zap("Example program to test conversion between Monte Carlo and DAQ times.");
85 
86  zap['f'] = make_field(inputFile);
87  zap['a'] = make_field(detectorFile);
88  zap['n'] = make_field(numberOfEvents) = JLimit::max();
89  zap['T'] = make_field(Tmax_ns) = 20.0;
90  zap['d'] = make_field(debug) = 3;
91 
92  zap(argc, argv);
93  }
94  catch(const exception &error) {
95  FATAL(error.what() << endl);
96  }
97 
98 
99  using namespace KM3NETDAQ;
100 
101  cout.tie(&cerr);
102 
103 
105 
106  try {
107  load(detectorFile, detector);
108  }
109  catch(const JException& error) {
110  FATAL(error);
111  }
112 
113  const JModuleRouter router(detector);
114 
115 
116  //typedef JDAQSnapshotHit JHit_t;
117  typedef JDAQTriggeredHit JHit_t;
118 
119  while (inputFile.hasNext()) {
120 
121  STATUS("event: " << setw(10) << inputFile.getCounter() << '\r'); DEBUG(endl);
122 
124 
125  const JDAQEvent* tev = ps;
126  const Evt* event = ps;
127 
128  ASSERT(!tev->empty<JHit_t>());
129 
130  const JTimeConverter converter(*event, *tev);
131 
132  typedef map<int, set<double> > map_type;
133 
134  map_type mc;
135  map_type data;
136 
137  for (vector<Hit>::const_iterator i = event->mc_hits.begin(); i != event->mc_hits.end(); ++i) {
138 
139  if (!is_noise(*i)) {
140 
141  DEBUG("Hit (MC): " << setw(5) << i->pmt_id << ' ' << FIXED(6,1) << getTime(*i) << endl);
142 
143  mc[i->pmt_id].insert(getTime(*i));
144  }
145  }
146 
147  for (JDAQEvent::const_iterator<JHit_t> i = tev->begin<JHit_t>(); i != tev->end<JHit_t>(); ++i) {
148 
149  const JModule& module = router.getModule(i->getModuleID());
150  const JPMT& pmt = module.getPMT(i->getPMT());
151 
152  DEBUG("Hit (DAQ): " << setw(5) << pmt.getID() << ' ' << FIXED(6,1) << converter.getTime(*i, pmt) << endl);
153 
154  data[pmt.getID()].insert(converter.getTime(*i, pmt));
155  }
156 
157  int number_of_hits = 0;
158 
159  for (map_type::const_iterator p = mc.begin(), q = data.begin(); ; ++p, ++q) {
160 
161  while (p != mc.end() && q != data.end() && p->first < q->first) { ++p; }
162  while (p != mc.end() && q != data.end() && q->first < p->first) { ++q; }
163 
164  if (p == mc.end() || q == data.end()) {
165  break;
166  }
167 
168  if (p->first == q->first) {
169  number_of_hits += get_count(p->second, q->second, Tmax_ns);
170  }
171  }
172 
173  NOTICE("Number of hits " << setw(5) << number_of_hits << '/' << setw(5) << tev->size<JHit_t>() << endl);
174 
175  ASSERT(number_of_hits != 0);
176  }
177 
178  ASSERT(inputFile.getCounter() != 0);
179 
180  return 0;
181 }
KM3NETDAQ::JDAQEvent
DAQ Event.
Definition: JDAQEvent.hh:30
JAANET::JHit_t
Auxiliary class to set-up Hit.
Definition: JHit_t.hh:25
JSUPPORT::JLimit
Auxiliary class for defining the range of iterations of objects.
Definition: JLimit.hh:41
FIXED
Auxiliary data structure for floating point format specification.
Definition: JPrint.hh:481
KM3NETDAQ::JDAQEvent::begin
const_iterator< T > begin() const
Get begin of data.
ASSERT
#define ASSERT(A,...)
Assert macro.
Definition: JMessage.hh:90
JDETECTOR::load
void load(const JString &file_name, JDetector &detector)
Load detector from input file.
Definition: JDetectorToolkit.hh:476
JLANG::JMultiPointer
General purpose class for multiple pointers.
Definition: JMultiPointer.hh:22
std::vector
Definition: JSTDTypes.hh:12
Evt
The Evt class respresent a Monte Carlo (MC) event as well as an offline event.
Definition: Evt.hh:19
JTRIGGER::JTimeConverter
Auxiliary class to convert DAQ/trigger hit time to/from Monte Carlo hit time.
Definition: JTimeConverter.hh:36
JPARSER::JParser
Utility class to parse command line options.
Definition: JParser.hh:1493
KM3NETDAQ::JDAQEvent::empty
bool empty() const
Check emptyness of hit container.
NOTICE
#define NOTICE(A)
Definition: JMessage.hh:64
KM3NETDAQ::JDAQEvent::end
const_iterator< T > end() const
Get end of data.
JSUPPORT::JTriggeredFileScanner::hasNext
virtual bool hasNext()
Check availability of next element.
Definition: JTriggeredFileScanner.hh:72
JPP
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
Definition: JAAnetToolkit.hh:37
KM3NETDAQ::JDAQTriggeredHit
DAQ triggered hit.
Definition: JDAQTriggeredHit.hh:20
KM3NETDAQ::JDAQEvent::size
unsigned int size() const
Get number of hits.
JSUPPORT::JMultipleFileScanner< JTypeList< JDAQEvent, JTypelist_t > >::getCounter
counter_type getCounter() const
Get counter.
Definition: JMultipleFileScanner.hh:323
JDETECTOR::JModule::getPMT
const JPMT & getPMT(const int index) const
Get PMT.
Definition: JModule.hh:173
debug
int debug
debug level
Definition: JSirene.cc:59
KM3NETDAQ::JDAQEvent::const_iterator
Template const_iterator.
Definition: JDAQEvent.hh:68
JSUPPORT::JTriggeredFileScanner::next
virtual const multi_pointer_type & next()
Get next element.
Definition: JTriggeredFileScanner.hh:92
JLANG::JObjectID::getID
int getID() const
Get identifier.
Definition: JObjectID.hh:55
JAANET::getTime
double getTime(const Hit &hit)
Get true time of hit.
Definition: JAAnetToolkit.hh:87
JDETECTOR::JModule
Data structure for a composite optical module.
Definition: JModule.hh:49
std::map
Definition: JSTDTypes.hh:16
JDETECTOR::JPMT
Data structure for PMT geometry and calibration.
Definition: JPMT.hh:53
STATUS
#define STATUS(A)
Definition: JMessage.hh:63
JDETECTOR::JModuleRouter
Router for direct addressing of module data in detector data structure.
Definition: JModuleRouter.hh:34
JSUPPORT::JTriggeredFileScanner
Auxiliary class to synchronously read DAQ events and Monte Carlo events (and optionally other events)...
Definition: JTriggeredFileScanner.hh:39
make_field
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1954
JDETECTOR::JDetector
Detector data structure.
Definition: JDetector.hh:80
JAANET::detector
Detector file.
Definition: JHead.hh:130
DEBUG
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62
std
Definition: jaanetDictionary.h:36
KM3NETDAQ
KM3NeT DAQ data structures and auxiliaries.
Definition: DataQueue.cc:39
FATAL
#define FATAL(A)
Definition: JMessage.hh:67
JLANG::JException
General exception.
Definition: JException.hh:23
JAANET::is_noise
bool is_noise(const Hit &hit)
Verify hit origin.
Definition: JAAnetToolkit.hh:121