Jpp 20.0.0-27-g39925593c-D
the software that should make you happy
Loading...
Searching...
No Matches
Functions
JTimeConverter.cc File Reference

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

#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 "km3net-dataformat/offline/io_online.hh"
#include "km3net-dataformat/tools/time_converter.hh"
#include "JDAQ/JDAQEventIO.hh"
#include "JAAnet/JAAnetToolkit.hh"
#include "JDetector/JDetector.hh"
#include "JDetector/JDetectorToolkit.hh"
#include "JDetector/JModuleRouter.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
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);
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
100
101 try {
103 }
104 catch(const JException& error) {
105 FATAL(error);
106 }
107
108 const JModuleRouter router(detector);
109
110
111 //typedef JDAQSnapshotHit JHit_t;
112 typedef JDAQTriggeredHit JHit_t;
113
114 while (inputFile.hasNext()) {
115
116 STATUS("event: " << setw(10) << inputFile.getCounter() << '\r'); DEBUG(endl);
117
119
120 const JDAQEvent* tev = ps;
121 const Evt* event = ps;
122
123 ASSERT(!tev->empty<JHit_t>());
124
126
127 map_type mc;
128 map_type data;
129
130 for (vector<Hit>::const_iterator i = event->mc_hits.begin(); i != event->mc_hits.end(); ++i) {
131
132 if (!is_noise(*i)) {
133
134 DEBUG("Hit (MC): " << setw(5) << i->pmt_id << ' ' << FIXED(6,1) << getTime(*i) << endl);
135
136 mc[i->pmt_id].insert(getTime(*i));
137 }
138 }
139
140 for (int test : {1, 2} ) {
141
142 DEBUG("Test "<< test << endl);
143
144 data.clear();
145
147
148 if (test == 1) {
149
150 converter = time_converter(*event, *tev);
151
152 for (JDAQEvent::const_iterator<JHit_t> i = tev->begin<JHit_t>(); i != tev->end<JHit_t>(); ++i) {
153
154 const JModule& module = router.getModule(i->getModuleID());
155 const JPMT& pmt = module.getPMT(i->getPMT());
156
157 data[pmt.getID()].insert(converter.getTime(getTime(*i, pmt)));
158 }
159
160 } else if (test == 2) {
161
162 Evt evt = *event;
163
164 read(evt, *tev);
165
167
168 for (vector<Hit>::const_iterator i = evt.hits.begin(); i != evt.hits.end(); ++i) {
169
170 if (i->trig != 0) {
171
172 const JModule& module = router.getModule(i->dom_id);
173 const JPMT& pmt = module.getPMT(i->channel_id);
174
175 data[pmt.getID()].insert(converter.getTime(getTime(i->tdc, pmt)));
176 }
177 }
178 }
179
180 for (map_type::const_iterator i = data.begin(); i != data.end(); ++i) {
181 for (set<double>::const_iterator hit = i->second.begin(); hit != i->second.end(); ++hit) {
182 DEBUG("Hit (DAQ): " << setw(5) << i->first << ' ' << FIXED(6,1) << *hit << endl);
183 }
184 }
185
186 int number_of_hits = 0;
187
188 for (map_type::const_iterator p = mc.begin(), q = data.begin(); ; ++p, ++q) {
189
190 while (p != mc.end() && q != data.end() && p->first < q->first) { ++p; }
191 while (p != mc.end() && q != data.end() && q->first < p->first) { ++q; }
192
193 if (p == mc.end() || q == data.end()) {
194 break;
195 }
196
197 if (p->first == q->first) {
198 number_of_hits += get_count(p->second, q->second, Tmax_ns);
199 }
200 }
201
202 NOTICE("Number of hits " << setw(5) << number_of_hits << '/' << setw(5) << tev->size<JHit_t>() << endl);
203
205 }
206 }
207
208 ASSERT(inputFile.getCounter() != 0);
209
210 return 0;
211}
#define DEBUG(A)
Message macros.
Definition JMessage.hh:62
#define STATUS(A)
Definition JMessage.hh:63
#define ASSERT(A,...)
Assert macro.
Definition JMessage.hh:90
#define NOTICE(A)
Definition JMessage.hh:64
#define FATAL(A)
Definition JMessage.hh:67
int debug
debug level
Definition JSirene.cc:72
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition JParser.hh:2142
Detector data structure.
Definition JDetector.hh:96
Router for direct addressing of module data in detector data structure.
Data structure for a composite optical module.
Definition JModule.hh:76
Data structure for PMT geometry, calibration and status.
Definition JPMT.hh:49
General exception.
Definition JException.hh:24
Template definition of a multi-dimensional oscillation probability interpolation table.
void insert(const JMultiFunction< __JFunction_t, __JMaplist_t, __JDistance_t > &input)
Insert multidimensional input.
Template const_iterator.
Definition JDAQEvent.hh:68
Auxiliary class to convert DAQ hit time to/from Monte Carlo hit time.
bool is_noise(const Hit &hit)
Verify hit origin.
std::istream & read(std::istream &in, JTestSummary &summary)
Read test summary.
void load(const std::string &file_name, JDetector &detector)
Load detector from input file.
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
const char * getTime()
Get current local time conform ISO-8601 standard.
KM3NeT DAQ data structures and auxiliaries.
Definition DataQueue.cc:39
std::map< int, range_type > map_type
The Evt class respresent a Monte Carlo (MC) event as well as an offline event.
Definition Evt.hh:21
std::vector< Hit > hits
list of hits
Definition Evt.hh:38
Auxiliary data structure for floating point format specification.
Definition JManip.hh:448
Detector file.
Definition JHead.hh:227
General purpose class for multiple pointers.
Auxiliary class to set-up Hit.
Definition JSirene.hh:60
Auxiliary class for defining the range of iterations of objects.
Definition JLimit.hh:45
static counter_type max()
Get maximum counter value.
Definition JLimit.hh:128