Jpp  18.0.0-rc.2
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Functions
JConvertDusj.cc File Reference
#include <string>
#include <iostream>
#include <iomanip>
#include <algorithm>
#include "TROOT.h"
#include "TFile.h"
#include "km3net-dataformat/offline/Head.hh"
#include "km3net-dataformat/offline/Evt.hh"
#include "km3net-dataformat/offline/io_online.hh"
#include "JDetector/JDetector.hh"
#include "JDetector/JDetectorToolkit.hh"
#include "JDetector/JModuleRouter.hh"
#include "JAAnet/JHead.hh"
#include "JAAnet/JHeadToolkit.hh"
#include "JDAQ/JDAQEventIO.hh"
#include "JDAQ/JDAQSummarysliceIO.hh"
#include "JSupport/JSingleFileScanner.hh"
#include "JSupport/JTreeScanner.hh"
#include "JSupport/JFileRecorder.hh"
#include "JSupport/JMonteCarloFileSupportkit.hh"
#include "JSupport/JSupportToolkit.hh"
#include "JSupport/JSupport.hh"
#include "JSupport/JMeta.hh"
#include "Jeep/JeepToolkit.hh"
#include "Jeep/JParser.hh"
#include "Jeep/JMessage.hh"

Go to the source code of this file.

Functions

int main (int argc, char **argv)
 Auxiliary program to convert data to Dusj format. More...
 

Function Documentation

int main ( int  argc,
char **  argv 
)

Auxiliary program to convert data to Dusj format.

Author
mdejong

Definition at line 44 of file JConvertDusj.cc.

45 {
46  using namespace std;
47  using namespace JPP;
48  using namespace KM3NETDAQ;
49 
52  JLimit_t& numberOfEvents = inputFile.getLimit();
53  string detectorFile;
54  int debug;
55 
56  try {
57 
58  JParser<> zap("Auxiliary program to convert data to Dusj format.");
59 
60  zap['f'] = make_field(inputFile);
61  zap['o'] = make_field(outputFile) = "dusj.evt";
62  zap['a'] = make_field(detectorFile);
63  zap['n'] = make_field(numberOfEvents) = JLimit::max();
64  zap['d'] = make_field(debug) = 2;
65 
66  zap(argc, argv);
67  }
68  catch(const exception& error) {
69  FATAL(error.what() << endl);
70  }
71 
72 
74 
75  try {
76  load(detectorFile, detector);
77  }
78  catch(const JException& error) {
79  FATAL(error);
80  }
81 
82  const JModuleRouter router(detector);
83 
84 
85  outputFile.open();
86 
87  Vec center(0,0,0);
88 
89  int mc_run_id = 0;
90 
91  Head header;
92 
93  try {
94 
95  header = getHeader(inputFile);
96 
97  JHead buffer(header);
98 
99  center = get<Vec>(buffer);
100  mc_run_id = buffer.start_run.run_id;
101 
102  } catch(const exception& error) {
103 
104  JHead buffer(header);
105 
106  buffer.start_run.run_id = getRunNumber<JDAQSummaryslice>(inputFile);
107  buffer.push(&JHead::start_run);
108 
109  copy(buffer, header);
110  }
111  {
112  JHead buffer(header);
113 
114  buffer.DAQ.livetime_s = getLivetime(inputFile);
115  buffer.push(&JHead::DAQ);
116 
117  copy(buffer, header);
118  }
119 
120  outputFile.put(header);
121  outputFile.put(JMeta(argc, argv));
122 
123  for (JTreeScanner<Evt> mc(inputFile); inputFile.hasNext(); ) {
124 
125  STATUS("event: " << setw(10) << inputFile.getCounter() << '\r'); DEBUG(endl);
126 
127  JDAQEvent* tev = inputFile.next();
128 
129  Evt out;
130 
131  if (mc.getEntries() != 0) {
132 
133  Evt* event = mc.getEntry(tev->getCounter());
134 
135  if (event != NULL) {
136 
137  out = *event;
138 
139  if (out.mc_run_id == 0) {
140  out.mc_run_id = mc_run_id;
141  }
142 
143  for (vector<Trk>::iterator i = out.mc_trks.begin(); i != out.mc_trks.end(); ++i) {
144  i->pos += center;
145  }
146  }
147 
148  } else {
149 
150  out.mc_id = inputFile.getCounter();
151  out.mc_run_id = tev->getRunNumber();
152  }
153 
154  read(out, *tev);
155 
156  out.hits.clear();
157 
159 
160  const JModule& module = router.getModule(i->getModuleID());
161  const JPMT& pmt = module.getPMT (i->getPMT());
162 
163  Hit hit;
164 
165  hit.id = out.hits.size() + 1;
166  hit.dom_id = i->getModuleID();
167  //hit.channel_id = i->getPMT();
168  //hit.tdc = i->getT();
169  //hit.tot = i->getToT();
170  hit.pmt_id = pmt.getID();
171  hit.t = i->getT(); // getTime(i->getT(), pmt.getCalibration());
172  hit.a = i->getToT(); // getToT (i->getToT(), pmt.getCalibration());
173  //hit.pos = Vec(pmt.getX(), pmt.getY(), pmt.getZ());
174  //hit.dir = Vec(pmt.getDX(), pmt.getDY(), pmt.getDZ());
175 
176  out.hits.push_back(hit);
177  }
178 
179  // time offset
180 
181  double t0 = numeric_limits<double>::max();
182 
183  for (vector<Hit>::const_iterator i = out.hits.begin(); i != out.hits.end(); ++i) {
184  if (t0 > i->t) {
185  t0 = i->t;
186  }
187  }
188 
189  for (vector<Hit>::iterator i = out.hits.begin(); i != out.hits.end(); ++i) {
190  i->t -= t0;
191  }
192 
193  if (out.mc_t > 0.0 && t0 != numeric_limits<double>::max()) {
194  out.mc_t -= t0;
195  }
196 
197  out.id = inputFile.getCounter();
198 
199  outputFile.put(out);
200  }
201 
202  STATUS(endl);
203  /*
204  JSingleFileScanner<JTYPELIST<JDAQSummaryslice, JMeta>::typelist> io(inputFile);
205 
206  io >> outputFile;
207  */
208  outputFile.close();
209 }
Auxiliary class for ROOT I/O of application specific meta data.
Definition: JMeta.hh:70
Object writing to file.
Utility class to parse command line options.
Definition: JParser.hh:1514
General exception.
Definition: JException.hh:23
DAQ key hit.
Definition: JDAQKeyHit.hh:19
Data structure for a composite optical module.
Definition: JModule.hh:68
std::istream & read(std::istream &in, JTestSummary &summary, const char delimiter= ' ')
Read test summary.
#define STATUS(A)
Definition: JMessage.hh:63
Detector data structure.
Definition: JDetector.hh:89
Template const_iterator.
Definition: JDAQEvent.hh:68
Router for direct addressing of module data in detector data structure.
int pmt_id
global PMT identifier as found in evt files
Definition: Hit.hh:20
string outputFile
int getRunNumber() const
Get run number.
Head getHeader(const JMultipleFileScanner_t &file_list)
Get Monte Carlo header.
Auxiliary class for defining the range of iterations of objects.
Definition: JLimit.hh:41
const_iterator< T > end() const
Get end of data.
double a
hit amplitude (in p.e.)
Definition: Hit.hh:24
The Vec class is a straightforward 3-d vector, which also works in pyroot.
Definition: Vec.hh:12
Detector file.
Definition: JHead.hh:226
int mc_run_id
MC run identifier.
Definition: Evt.hh:27
const_iterator< T > begin() const
Get begin of data.
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1989
int mc_id
identifier of the MC event (as found in ascii or antcc file).
Definition: Evt.hh:24
int getID() const
Get identifier.
Definition: JObjectID.hh:50
Data structure for PMT geometry, calibration and status.
Definition: JPMT.hh:43
Q DAQ
Definition: JDataQuality.sh:59
int id
Definition: Hit.hh:11
const JPMT & getPMT(const int index) const
Get PMT.
Definition: JModule.hh:173
double mc_t
MC: time where the mc-event was put in the timeslice, since start of run (offset+frameidx*timeslice_d...
Definition: Evt.hh:47
Monte Carlo run header.
Definition: JHead.hh:1221
The Head class reflects the header of Monte-Carlo event files, which consists of keys (also referred ...
Definition: Head.hh:65
#define FATAL(A)
Definition: JMessage.hh:67
Definition: Hit.hh:8
void load(const std::string &file_name, JDetector &detector)
Load detector from input file.
double getLivetime(const std::string &file_name)
Get data taking live time.
double t
hit time (from tdc+calibration or MC truth)
Definition: Hit.hh:23
int dom_id
module identifier from the data (unique in the detector).
Definition: Hit.hh:14
void copy(const Head &from, JHead &to)
Copy header from from to to.
Definition: JHead.cc:162
int id
offline event identifier
Definition: Evt.hh:22
Object reading from a list of files.
const JLimit & getLimit() const
Get limit.
Definition: JLimit.hh:84
std::vector< Hit > hits
list of hits
Definition: Evt.hh:38
do set_variable DETECTOR_TXT $WORKDIR detector
JTriggerCounter_t getCounter() const
Get trigger counter.
std::vector< Trk > mc_trks
MC: list of MC truth tracks.
Definition: Evt.hh:49
int debug
debug level
JTriggerCounter_t next()
Increment trigger counter.
The Evt class respresent a Monte Carlo (MC) event as well as an offline event.
Definition: Evt.hh:20
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62