Jpp 20.0.0-195-g190c9e876
the software that should make you happy
Loading...
Searching...
No Matches
software/JReconstruction/JConvertEvt.cc
Go to the documentation of this file.
1
2#include <string>
3#include <iostream>
4#include <iomanip>
5#include <algorithm>
6#include <memory>
7#include <set>
8#include <map>
9
10#include "TROOT.h"
11#include "TFile.h"
12
18
25
27
28#include "JCompass/JEvt.hh"
29#include "JCompass/JSupport.hh"
30
31#include "JAcoustics/JEvt.hh"
33
35
36#include "JAAnet/JHead.hh"
39
40#include "JDAQ/JDAQEventIO.hh"
42
48#include "JSupport/JSupport.hh"
49#include "JSupport/JMeta.hh"
50
54
55#include "Jeep/JeepToolkit.hh"
56#include "Jeep/JParser.hh"
57#include "Jeep/JMessage.hh"
58
59
60/**
61 * \file
62 * Auxiliary program to convert fit results to Evt format.
63 *
64 * \author mdejong
65 */
66int main(int argc, char **argv)
67{
68 using namespace std;
69 using namespace JPP;
70 using namespace KM3NETDAQ;
71
72 typedef JMultipleFileScanner<calibration_types> JCalibration_t;
73
74 JMultipleFileScanner<> inputFile;
76 JLimit_t& numberOfEvents = inputFile.getLimit();
77 string detectorFile;
78 JCalibration_t calibrationFile;
79 double Tmax_s;
80 JPMTParametersMap pmtParameters;
82 size_t numberOfFits;
83 int debug;
84
85 try {
86
87 JParser<> zap("Auxiliary program to convert fit results to Evt format.\
88 \nThe option -L corresponds to the name of a shared library \
89 \nand function so to rearrange the order of fit results.");
90
91 zap['f'] = make_field(inputFile);
92 zap['o'] = make_field(outputFile) = "aanet.root";
93 zap['n'] = make_field(numberOfEvents) = JLimit::max();
94 zap['a'] = make_field(detectorFile) = "";
95 zap['+'] = make_field(calibrationFile) = JPARSER::initialised();
96 zap['T'] = make_field(Tmax_s) = 100.0;
97 zap['P'] = make_field(pmtParameters, "PMT simulation data (or corresponding file name)") = JPARSER::initialised();
99 zap['N'] = make_field(numberOfFits) = 0;
100 zap['d'] = make_field(debug) = 2;
101
102 zap(argc, argv);
103 }
104 catch(const exception& error) {
105 FATAL(error.what() << endl);
106 }
107
108 if (detectorFile == "" && !calibrationFile.empty()) {
109 FATAL("Missing detector file." << endl);
110 }
111
112 typedef JParallelFileScanner< JTypeList<JDAQEvent, JFIT::JEvt> > JParallelFileScanner_t;
113 typedef JParallelFileScanner_t::multi_pointer_type multi_pointer_type;
114
116
117 if (detectorFile != "") {
118 try {
119 load(detectorFile, detector);
120 }
121 catch(const JException& error) {
122 FATAL(error);
123 }
124 }
125
126 unique_ptr<JDynamics> dynamics;
127
128 if (!calibrationFile.empty()) {
129
130 try {
131
132 dynamics.reset(new JDynamics(detector, Tmax_s));
133
134 dynamics->load(calibrationFile);
135 }
136 catch(const exception& error) {
137 FATAL(error.what());
138 }
139 }
140
141
142 const JPMTRouter pmt_router(dynamics ? dynamics->getDetector() : detector);
143 const JModuleRouter mod_router(dynamics ? dynamics->getDetector() : detector);
144
145
146 outputFile.open();
147
148 Head header;
149
150 try {
151 header = getHeader(inputFile);
152 } catch(const exception& error) {}
153
154 JAANET::JHead buffer(header);
155
156 // DAQ livetime is not yet set for real data, for Monte Carlo, it is set in JTriggerEfficiency
157
158 if (buffer.pull(&JAANET::JHead::DAQ) == buffer.end()) {
159
160 buffer.DAQ.livetime_s = getLivetime(inputFile.begin(), inputFile.end());
161 buffer.push(&JAANET::JHead::DAQ);
162
163 copy(buffer, header);
164 }
165
166 if (detectorFile != "") {
167
170
171 copy(buffer, header);
172 }
173
174
175 outputFile.put(header);
176
177
178 outputFile.put(JMeta(argc, argv));
179
180 map<int, size_t> miss_mod;
181 map<int, size_t> miss_pmt;
182
183 counter_type number_of_events = 0;
184
185 for (JMultipleFileScanner<>::const_iterator i = inputFile.begin(); i != inputFile.end(); ++i) {
186
187 STATUS("Processing: " << *i << endl);
188
189 JParallelFileScanner_t in(*i);
190 JTreeScanner<Evt> mc(*i);
191
192 in.setLimit(inputFile.getLimit());
193
194 Vec offset(0,0,0);
195
196 int mc_run_id = 0;
197
198 try {
199
200 const JAANET::JHead head(getHeader(*i));
201
202 offset = getOffset(head);
203 mc_run_id = head.start_run.run_id;
204
205 } catch(const exception& error) {}
206
207
208 while (in.hasNext()) {
209
210 STATUS("event: " << setw(10) << in.getCounter() << '\r'); DEBUG(endl);
211
212 multi_pointer_type ps = in.next();
213
214 JDAQEvent* tev = ps;
215 JFIT::JEvt* evt = ps;
216
217 if (dynamics) {
218 dynamics->update(*tev); // update detector
219 }
220
221 JFIT::JEvt::iterator __end = evt->end();
222
223 if (numberOfFits > 0) {
224 advance(__end = evt->begin(), min(numberOfFits, evt->size()));
225 }
226
227 if (qualitySorter.is_valid()) {
228 partial_sort(evt->begin(), __end, evt->end(), qualitySorter);
229 }
230
231 Evt out; // output event
232
233 if (mc.getEntries() != 0) {
234
235 Evt* event = mc.getEntry(tev->getCounter());
236
237 if (event != NULL) {
238
239 out = *event; // import Monte Carlo true data
240
241 if (out.mc_run_id == 0) {
242 out.mc_run_id = mc_run_id;
243 }
244
245 for (vector<Trk>::iterator i = out.mc_trks.begin(); i != out.mc_trks.end(); ++i) {
246 i->pos += offset; // offset true positions
247 }
248 }
249 }
250
251 read(out, *tev); // import DAQ data
252
253 if (!pmt_router->empty()) { // import calibration data
254
255 for (vector<Hit>::iterator i = out.mc_hits.begin(); i != out.mc_hits.end(); ++i) {
256
257 if (pmt_router.hasPMT(i->pmt_id)) {
258
259 const JPMTAddress address = pmt_router.getAddress(i->pmt_id);
260 const JPMTIdentifier id = pmt_router.getIdentifier(address);
261 const JPMT& pmt = pmt_router.getPMT(address);
262
263 i->dom_id = id.getID();
264 i->channel_id = id.getPMTAddress();
265 i->pos = getPosition (pmt);
266 i->dir = getDirection(pmt);
267
268 } else {
269
270 miss_pmt[i->pmt_id] += 1;
271 }
272 }
273 }
274
275 if (!mod_router->empty()) { // import calibration data
276
277 for (vector<Hit>::iterator i = out.hits.begin(); i != out.hits.end(); ++i) {
278
279 if (mod_router.hasModule(i->dom_id)) {
280
281 const JPMTIdentifier id(i->dom_id, i->channel_id);
282
283 const JPMTParameters& parameters = pmtParameters.getPMTParameters(id);
284 const JPMTAnalogueSignalProcessor cpu(parameters);
285
286 const JPMT& pmt = mod_router.getPMT(id);
287 const int type = parameters.getType();
288
289 i->pmt_id = pmt.getID();
290 i->pos = getPosition (pmt);
291 i->dir = getDirection(pmt);
292 i->t = getTime(i->tdc, pmt.getCalibration()) - (parameters.slewing ? getTimeSlewing(i->tot, type) : 0.0);
293 i->a = cpu.getNPE(i->tot);
294
295 } else {
296
297 miss_mod[i->dom_id] += 1;
298 }
299 }
300 }
301
302 copy(evt->cbegin(), __end, out);
303
304 out.id = ++number_of_events;
305
306 outputFile.put(out);
307 }
308 }
309
310 STATUS(endl);
311
312 for (const auto& i : miss_pmt) { ERROR("Misses PMT " << setw(8) << i.first << ' ' << setw(8) << i.second << endl); }
313 for (const auto& i : miss_mod) { ERROR("Misses module " << setw(8) << i.first << ' ' << setw(8) << i.second << endl); }
314
316
317 io >> outputFile;
318
319 outputFile.close();
320}
Definition of hit and track types and auxiliary methods for handling Monte Carlo data.
Acoustic event fit.
ROOT TTree parameter settings.
Compass event data types.
ROOT TTree parameter settings.
string outputFile
Data structure for detector geometry and calibration.
Dynamic detector calibration.
Recording of objects on file according a format that follows from the file name extension.
General purpose messaging.
#define DEBUG(A)
Message macros.
Definition JMessage.hh:62
#define STATUS(A)
Definition JMessage.hh:63
#define FATAL(A)
Definition JMessage.hh:67
int debug
debug level
Definition JSirene.cc:74
ROOT I/O of application specific meta data.
Direct access to module in detector data structure.
PMT analogue signal processor.
Direct access to PMT in detector data structure.
Parallel scanning of objects from a single file or multiple files according a format that follows fro...
Utility class to parse command line options.
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition JParser.hh:2142
Support methods.
ROOT TTree parameter settings of various packages.
Time slewing correction.
Auxiliary methods for handling file names, type names and environment.
Monte Carlo run header.
Definition JHead.hh:1236
JAANET::start_run start_run
Definition JHead.hh:1601
JAANET::DAQ DAQ
Definition JHead.hh:1625
void push(T JHead::*pd)
Push given data member to Head.
Definition JHead.hh:1374
const_iterator pull(T JHead::*pd) const
Pull given data member from Head.
Definition JHead.hh:1349
JAANET::calibration calibration
Definition JHead.hh:1610
const JCalibration & getCalibration() const
Get calibration.
Detector data structure.
Definition JDetector.hh:96
Router for direct addressing of module data in detector data structure.
const JPMT & getPMT(const JPMTIdentifier &id) const
Get PMT parameters.
bool hasModule(const JObjectID &id) const
Has module.
Address of PMT in detector data structure.
Auxiliary class for map of PMT parameters.
const JPMTParameters & getPMTParameters(const JPMTIdentifier &id) const
Get PMT parameters.
Data structure for PMT parameters.
int getType() const
Get type for for time-slewing correction.
bool slewing
time slewing of analogue signal
Router for direct addressing of PMT data in detector data structure.
Definition JPMTRouter.hh:37
bool hasPMT(const JObjectID &id) const
Has PMT.
const JPMT & getPMT(const JPMTAddress &address) const
Get PMT.
Definition JPMTRouter.hh:92
JPMTIdentifier getIdentifier(const JPMTAddress &address) const
Get identifier of PMT.
const JPMTAddress & getAddress(const JObjectID &id) const
Get address of PMT.
Definition JPMTRouter.hh:80
Data structure for PMT geometry, calibration and status.
Definition JPMT.hh:49
Data structure for set of track fit results.
General exception.
Definition JException.hh:24
int getID() const
Get identifier.
Definition JObjectID.hh:50
Utility class to parse command line options.
Definition JParser.hh:1698
Object writing to file.
General purpose class for object reading from a list of file names.
General purpose class for parallel reading of objects from a single file or multiple files.
Template definition for direct access of elements in ROOT TChain.
JTriggerCounter_t getCounter() const
Get trigger counter.
JDirection3D getDirection(const Vec &dir)
Get direction.
void copy(const Head &from, JHead &to)
Copy header from from to to.
Definition JHead.cc:163
JPosition3D getPosition(const Vec &pos)
Get position.
Vec getOffset(const JHead &header)
Get offset.
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).
bool qualitySorter(const JFit &first, const JFit &second)
Comparison of fit results.
counter_type advance(counter_type &counter, const counter_type value, const counter_type limit=std::numeric_limits< counter_type >::max())
Advance counter.
Long64_t counter_type
Type definition for counter.
double getLivetime(const std::string &file_name)
Get data taking live time.
Head getHeader(const JMultipleFileScanner_t &file_list)
Get Monte Carlo header.
const char * getTime()
Get current local time conform ISO-8601 standard.
JTRIGGER::JTimeSlewing_t getTimeSlewing
Function object to get time-slewing correction.
KM3NeT DAQ data structures and auxiliaries.
Definition DataQueue.cc:39
int main(int argc, char **argv)
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
int mc_run_id
MC run identifier.
Definition Evt.hh:27
std::vector< Hit > mc_hits
MC: list of MC truth hits.
Definition Evt.hh:48
std::vector< Trk > mc_trks
MC: list of MC truth tracks.
Definition Evt.hh:49
int id
offline event identifier
Definition Evt.hh:22
The Head class reflects the header of Monte-Carlo event files, which consists of keys (also referred ...
Definition Head.hh:65
double livetime_s
Live time [s].
Definition JHead.hh:1053
std::string buffer
General purpose name.
Definition JHead.hh:217
static const std::string dynamical()
Definition JHead.hh:332
static const std::string statical()
Definition JHead.hh:331
Detector file.
Definition JHead.hh:227
int run_id
MC run number.
Definition JHead.hh:143
virtual double getNPE(const double tot_ns) const override
Get number of photo-electrons.
Dynamic detector calibration.
Definition JDynamics.hh:81
Empty structure for specification of parser element that is initialised (i.e. does not require input)...
Definition JParser.hh:68
General purpose sorter of fit results.
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
Auxiliary class for ROOT I/O of application specific meta data.
Definition JMeta.hh:72
The Vec class is a straightforward 3-d vector, which also works in pyroot.
Definition Vec.hh:13