Jpp 20.0.0-rc.2
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
34#include "JTrigger/JHit.hh"
36
37#include "JAAnet/JHead.hh"
40
41#include "JDAQ/JDAQEventIO.hh"
43
49#include "JSupport/JSupport.hh"
50#include "JSupport/JMeta.hh"
51
55
56#include "Jeep/JeepToolkit.hh"
57#include "Jeep/JParser.hh"
58#include "Jeep/JMessage.hh"
59
60
61/**
62 * \file
63 * Auxiliary program to convert fit results to Evt format.
64 *
65 * \author mdejong
66 */
67int main(int argc, char **argv)
68{
69 using namespace std;
70 using namespace JPP;
71 using namespace KM3NETDAQ;
72
74 JACOUSTICS::JEvt>::typelist calibration_types;
75 typedef JMultipleFileScanner<calibration_types> JCalibration_t;
76
77 JMultipleFileScanner<> inputFile;
79 JLimit_t& numberOfEvents = inputFile.getLimit();
80 string detectorFile;
81 JCalibration_t calibrationFile;
82 double Tmax_s;
83 JPMTParametersMap pmtParameters;
85 size_t numberOfFits;
86 int debug;
87
88 try {
89
90 JParser<> zap("Auxiliary program to convert fit results to Evt format.\
91 \nThe option -L corresponds to the name of a shared library \
92 \nand function so to rearrange the order of fit results.");
93
94 zap['f'] = make_field(inputFile);
95 zap['o'] = make_field(outputFile) = "aanet.root";
96 zap['n'] = make_field(numberOfEvents) = JLimit::max();
97 zap['a'] = make_field(detectorFile) = "";
98 zap['+'] = make_field(calibrationFile) = JPARSER::initialised();
99 zap['T'] = make_field(Tmax_s) = 100.0;
100 zap['P'] = make_field(pmtParameters, "PMT simulation data (or corresponding file name)") = JPARSER::initialised();
102 zap['N'] = make_field(numberOfFits) = 0;
103 zap['d'] = make_field(debug) = 2;
104
105 zap(argc, argv);
106 }
107 catch(const exception& error) {
108 FATAL(error.what() << endl);
109 }
110
111 if (detectorFile == "" && !calibrationFile.empty()) {
112 FATAL("Missing detector file." << endl);
113 }
114
115 typedef JParallelFileScanner< JTypeList<JDAQEvent, JFIT::JEvt> > JParallelFileScanner_t;
116 typedef JParallelFileScanner_t::multi_pointer_type multi_pointer_type;
117
119
120 if (detectorFile != "") {
121 try {
122 load(detectorFile, detector);
123 }
124 catch(const JException& error) {
125 FATAL(error);
126 }
127 }
128
129 unique_ptr<JDynamics> dynamics;
130
131 if (!calibrationFile.empty()) {
132
133 try {
134
135 dynamics.reset(new JDynamics(detector, Tmax_s));
136
137 dynamics->load(calibrationFile);
138 }
139 catch(const exception& error) {
140 FATAL(error.what());
141 }
142 }
143
144
145 const JPMTRouter pmt_router(dynamics ? dynamics->getDetector() : detector);
146 const JModuleRouter mod_router(dynamics ? dynamics->getDetector() : detector);
147
148
149 outputFile.open();
150
151 Head header;
152
153 try {
154 header = getHeader(inputFile);
155 } catch(const exception& error) {}
156
157 JAANET::JHead buffer(header);
158
159 // DAQ livetime is not yet set for real data, for Monte Carlo, it is set in JTriggerEfficiency
160
161 if (buffer.pull(&JAANET::JHead::DAQ) == buffer.end()) {
162
163 buffer.DAQ.livetime_s = getLivetime(inputFile.begin(), inputFile.end());
164 buffer.push(&JAANET::JHead::DAQ);
165
166 copy(buffer, header);
167 }
168
169 if (detectorFile != "") {
170
173
174 copy(buffer, header);
175 }
176
177
178 outputFile.put(header);
179
180
181 outputFile.put(JMeta(argc, argv));
182
183 map<int, size_t> miss_mod;
184 map<int, size_t> miss_pmt;
185
186 counter_type number_of_events = 0;
187
188 for (JMultipleFileScanner<>::const_iterator i = inputFile.begin(); i != inputFile.end(); ++i) {
189
190 STATUS("Processing: " << *i << endl);
191
192 JParallelFileScanner_t in(*i);
193 JTreeScanner<Evt> mc(*i);
194
195 in.setLimit(inputFile.getLimit());
196
197 Vec offset(0,0,0);
198
199 int mc_run_id = 0;
200
201 try {
202
203 const JAANET::JHead head(getHeader(*i));
204
205 offset = getOffset(head);
206 mc_run_id = head.start_run.run_id;
207
208 } catch(const exception& error) {}
209
210
211 while (in.hasNext()) {
212
213 STATUS("event: " << setw(10) << in.getCounter() << '\r'); DEBUG(endl);
214
215 multi_pointer_type ps = in.next();
216
217 JDAQEvent* tev = ps;
218 JFIT::JEvt* evt = ps;
219
220 if (dynamics) {
221 dynamics->update(*tev); // update detector
222 }
223
224 JFIT::JEvt::iterator __end = evt->end();
225
226 if (numberOfFits > 0) {
227 advance(__end = evt->begin(), min(numberOfFits, evt->size()));
228 }
229
230 if (qualitySorter.is_valid()) {
231 partial_sort(evt->begin(), __end, evt->end(), qualitySorter);
232 }
233
234 Evt out; // output event
235
236 if (mc.getEntries() != 0) {
237
238 Evt* event = mc.getEntry(tev->getCounter());
239
240 if (event != NULL) {
241
242 out = *event; // import Monte Carlo true data
243
244 if (out.mc_run_id == 0) {
245 out.mc_run_id = mc_run_id;
246 }
247
248 for (vector<Trk>::iterator i = out.mc_trks.begin(); i != out.mc_trks.end(); ++i) {
249 i->pos += offset; // offset true positions
250 }
251 }
252 }
253
254 read(out, *tev); // import DAQ data
255
256 if (!pmt_router->empty()) { // import calibration data
257
258 for (vector<Hit>::iterator i = out.mc_hits.begin(); i != out.mc_hits.end(); ++i) {
259
260 if (pmt_router.hasPMT(i->pmt_id)) {
261
262 const JPMTAddress address = pmt_router.getAddress(i->pmt_id);
263 const JPMTIdentifier id = pmt_router.getIdentifier(address);
264 const JPMT& pmt = pmt_router.getPMT(address);
265
266 i->dom_id = id.getID();
267 i->channel_id = id.getPMTAddress();
268 i->pos = getPosition (pmt);
269 i->dir = getDirection(pmt);
270
271 } else {
272
273 miss_pmt[i->pmt_id] += 1;
274 }
275 }
276 }
277
278 if (!mod_router->empty()) { // import calibration data
279
280 for (vector<Hit>::iterator i = out.hits.begin(); i != out.hits.end(); ++i) {
281
282 if (mod_router.hasModule(i->dom_id)) {
283
284 const JPMTIdentifier id(i->dom_id, i->channel_id);
285
286 const JPMTParameters& parameters = pmtParameters.getPMTParameters(id);
287 const JPMTAnalogueSignalProcessor cpu(parameters);
288
289 const JPMT& pmt = mod_router.getPMT(id);
290
292
293 const JTRIGGER::JHit hit(getTime(i->tdc, pmt.getCalibration()), i->tot);
294
295 i->pmt_id = pmt.getID();
296 i->pos = getPosition (pmt);
297 i->dir = getDirection(pmt);
298 i->t = hit.getT();
299 i->tot = hit.getToT();
300 i->a = cpu.getNPE(i->tot);
301
302 } else {
303
304 miss_mod[i->dom_id] += 1;
305 }
306 }
307 }
308
309 copy(evt->cbegin(), __end, out);
310
311 out.id = ++number_of_events;
312
313 outputFile.put(out);
314 }
315 }
316
317 STATUS(endl);
318
319 for (const auto& i : miss_pmt) { ERROR("Misses PMT " << setw(8) << i.first << ' ' << setw(8) << i.second << endl); }
320 for (const auto& i : miss_mod) { ERROR("Misses module " << setw(8) << i.first << ' ' << setw(8) << i.second << endl); }
321
323
324 io >> outputFile;
325
326 outputFile.close();
327}
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.
Tools for handling different hit types.
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:72
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.
Basic data structure for time and time over threshold information of hit.
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.
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.
Hit data structure.
static void setSlewing(const bool slewing)
Set slewing option.
double getToT() const
Get calibrated time over threshold of hit.
double getT() const
Get calibrated time of hit.
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.
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
Acoustic event fit.
Orientation of module.
virtual double getNPE(const double tot_ns) const override
Get number of photo-electrons.
Dynamic detector calibration.
Definition JDynamics.hh:86
Auxiliary class for recursive type list generation.
Definition JTypeList.hh:351
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