Jpp  15.0.0
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Functions
software/JReconstruction/JConvertEvt.cc File Reference

Auxiliary program to convert fit results to Evt format. More...

#include <string>
#include <iostream>
#include <iomanip>
#include <algorithm>
#include <memory>
#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/JPMTRouter.hh"
#include "JDetector/JModuleRouter.hh"
#include "JDetector/JPMTParametersMap.hh"
#include "JDetector/JPMTAnalogueSignalProcessor.hh"
#include "JDynamics/JDynamics.hh"
#include "JCompass/JEvt.hh"
#include "JCompass/JSupport.hh"
#include "JAcoustics/JEvt.hh"
#include "JAcoustics/JSupport.hh"
#include "JTrigger/JHit.hh"
#include "JTrigger/JHitToolkit.hh"
#include "JAAnet/JHead.hh"
#include "JAAnet/JHeadToolkit.hh"
#include "JAAnet/JAAnetToolkit.hh"
#include "JDAQ/JDAQEventIO.hh"
#include "JDAQ/JDAQSummarysliceIO.hh"
#include "JSupport/JParallelFileScanner.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 "JReconstruction/JEvt.hh"
#include "JReconstruction/JEvtToolkit.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)
 

Detailed Description

Auxiliary program to convert fit results to Evt format.

Author
mdejong

Definition in file software/JReconstruction/JConvertEvt.cc.

Function Documentation

int main ( int  argc,
char **  argv 
)

Definition at line 62 of file software/JReconstruction/JConvertEvt.cc.

63 {
64  using namespace std;
65  using namespace JPP;
66  using namespace KM3NETDAQ;
67 
69  JACOUSTICS::JEvt>::typelist calibration_types;
70  typedef JMultipleFileScanner<calibration_types> JCalibration_t;
71 
72  JMultipleFileScanner<> inputFile;
74  JLimit_t& numberOfEvents = inputFile.getLimit();
75  string detectorFile;
76  JCalibration_t calibrationFile;
77  double Tmax_s;
78  JPMTParametersMap pmtParameters;
80  size_t numberOfFits;
81  int debug;
82 
83  try {
84 
85  JParser<> zap("Auxiliary program to convert fit results to Evt format.\
86  \nThe option -L corresponds to the name of a shared library \
87  \nand function so to rearrange the order of fit results.");
88 
89  zap['f'] = make_field(inputFile);
90  zap['o'] = make_field(outputFile) = "aanet.root";
91  zap['n'] = make_field(numberOfEvents) = JLimit::max();
92  zap['a'] = make_field(detectorFile) = "";
93  zap['+'] = make_field(calibrationFile) = JPARSER::initialised();
94  zap['T'] = make_field(Tmax_s) = 100.0;
95  zap['P'] = make_field(pmtParameters, "PMT simulation data (or corresponding file name)") = JPARSER::initialised();
97  zap['N'] = make_field(numberOfFits) = 0;
98  zap['d'] = make_field(debug) = 2;
99 
100  zap(argc, argv);
101  }
102  catch(const exception& error) {
103  FATAL(error.what() << endl);
104  }
105 
106 
107  typedef JParallelFileScanner< JTypeList<JDAQEvent, JFIT::JEvt> > JParallelFileScanner_t;
108  typedef JParallelFileScanner_t::multi_pointer_type multi_pointer_type;
109 
110 
112 
113  if (detectorFile != "") {
114  try {
115  load(detectorFile, detector);
116  }
117  catch(const JException& error) {
118  FATAL(error);
119  }
120  }
121 
122  getMechanics.load(detector.getID());
123 
124  unique_ptr<JDynamics> dynamics;
125 
126  try {
127 
128  dynamics.reset(new JDynamics(detector, Tmax_s));
129 
130  dynamics->load(calibrationFile);
131  }
132  catch(const exception& error) {
133  if (!calibrationFile.empty()) {
134  FATAL(error.what());
135  }
136  }
137 
138 
139  const JPMTRouter pmt_router(dynamics ? dynamics->getDetector() : detector);
140  const JModuleRouter mod_router(dynamics ? dynamics->getDetector() : detector);
141 
142 
143  outputFile.open();
144 
145  Head header;
146 
147  try {
148  header = getHeader(inputFile);
149  } catch(const exception& error) {}
150  {
151  JAANET::JHead buffer(header);
152 
153  buffer.DAQ.livetime_s = getLivetime(inputFile.begin(), inputFile.end());
154  buffer.push(&JAANET::JHead::DAQ);
155 
156  copy(buffer, header);
157  }
158 
159  outputFile.put(header);
160 
161 
162  outputFile.put(JMeta(argc, argv));
163 
164 
165  counter_type number_of_events = 0;
166 
167  for (JMultipleFileScanner<>::const_iterator i = inputFile.begin(); i != inputFile.end(); ++i) {
168 
169  STATUS("Processing: " << *i << endl);
170 
171  JParallelFileScanner_t in(*i);
172  JTreeScanner<Evt> mc(*i);
173 
174  in.setLimit(inputFile.getLimit());
175 
176  Vec center(0,0,0);
177 
178  int mc_run_id = 0;
179 
180  try {
181 
182  const JAANET::JHead head(getHeader(*i));
183 
184  center = get<Vec>(head);
185  mc_run_id = head.start_run.run_id;
186 
187  } catch(const exception& error) {}
188 
189 
190  while (in.hasNext()) {
191 
192  STATUS("event: " << setw(10) << in.getCounter() << '\r'); DEBUG(endl);
193 
194  multi_pointer_type ps = in.next();
195 
196  JDAQEvent* tev = ps;
197  JFIT::JEvt* evt = ps;
198 
199  if (dynamics) {
200  dynamics->update(*tev); // update detector
201  }
202 
203  JFIT::JEvt::iterator __end = evt->end();
204 
205  if (numberOfFits > 0) {
206  advance(__end = evt->begin(), min(numberOfFits, evt->size()));
207  }
208 
209  if (qualitySorter.is_valid()) {
210  partial_sort(evt->begin(), __end, evt->end(), qualitySorter);
211  }
212 
213  Evt out; // output event
214 
215  if (mc.getEntries() != 0) {
216 
217  Evt* event = mc.getEntry(tev->getCounter());
218 
219  if (event != NULL) {
220 
221  out = *event; // import Monte Carlo true data
222 
223  if (out.mc_run_id == 0) {
224  out.mc_run_id = mc_run_id;
225  }
226 
227  for (vector<Trk>::iterator i = out.mc_trks.begin(); i != out.mc_trks.end(); ++i) {
228  i->pos += center; // offset true positions
229  }
230  }
231  }
232 
233  read(out, *tev); // import DAQ data
234 
235  if (!pmt_router->empty()) { // import calibration data
236 
237  for (vector<Hit>::iterator i = out.mc_hits.begin(); i != out.mc_hits.end(); ++i) {
238 
239  if (pmt_router.hasPMT(i->pmt_id)) {
240 
241  const JPMTAddress address = pmt_router.getAddress(i->pmt_id);
242  const JPMTIdentifier id = pmt_router.getIdentifier(address);
243  const JPMT& pmt = pmt_router.getPMT(address);
244 
245  i->dom_id = id.getID();
246  i->channel_id = id.getPMTAddress();
247  i->pos = getPosition (pmt);
248  i->dir = getDirection(pmt);
249 
250  } else {
251 
252  FATAL("Missing PMT" << i->pmt_id << endl);
253  }
254  }
255  }
256 
257  if (!mod_router->empty()) { // import calibration data
258 
259  for (vector<Hit>::iterator i = out.hits.begin(); i != out.hits.end(); ++i) {
260 
261  if (mod_router.hasModule(i->dom_id)) {
262 
263  const JPMTIdentifier id(i->dom_id, i->channel_id);
264 
265  const JPMTParameters& parameters = pmtParameters.getPMTParameters(id);
266  const JPMTAnalogueSignalProcessor cpu(parameters);
267 
268  const JPMT& pmt = mod_router.getPMT(id);
269 
271 
272  const JTRIGGER::JHit hit(getTime(i->tdc, pmt.getCalibration()), i->tot);
273 
274  i->pmt_id = pmt.getID();
275  i->pos = getPosition (pmt);
276  i->dir = getDirection(pmt);
277  i->t = hit.getT();
278  i->tot = hit.getToT();
279  i->a = cpu.getNPE(i->tot);
280 
281  } else {
282 
283  FATAL("Missing module " << i->dom_id << endl);
284  }
285  }
286  }
287 
288  for (JFIT::JEvt::const_iterator fit = evt->begin(); fit != __end; ++fit) { // import reconstruction data
289 
290  Trk trk;
291 
292  trk.id = out.trks.size() + 1;
293  trk.pos = Vec(fit->getX(), fit->getY(), fit->getZ());
294  trk.dir = Vec(fit->getDX(), fit->getDY(), fit->getDZ());
295  trk.t = fit->getT();
296  trk.E = fit->getE();
297  trk.lik = fit->getQ();
299 
300  for (JHistory::const_iterator i = fit->getHistory().begin(); i != fit->getHistory().end(); ++i) {
301  trk.rec_stages.push_back(i->type);
302  }
303 
304  for (int i = 0; i != fit->getN(); ++i) {
305  trk.fitinf.push_back(fit->getW(i));
306  }
307 
308  trk.error_matrix = fit->getV();
309 
310  out.trks.push_back(trk);
311  }
312 
313  out.id = ++number_of_events;
314 
315  outputFile.put(out);
316  }
317  }
318 
319  STATUS(endl);
320 
322 
323  io >> outputFile;
324 
325  outputFile.close();
326 }
Router for direct addressing of PMT data in detector data structure.
Definition: JPMTRouter.hh:33
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:1500
General exception.
Definition: JException.hh:23
std::istream & read(std::istream &in, JTestSummary &summary, const char delimiter= ' ')
Read test summary.
static void setSlewing(const bool slewing)
Set slewing option.
double t
track time [ns] (when the particle is at pos )
Definition: Trk.hh:19
static JDetectorMechanics getMechanics
Function object to get string mechanics.
Definition: JMechanics.hh:243
const JCalibration & getCalibration() const
Get calibration.
#define STATUS(A)
Definition: JMessage.hh:63
Detector data structure.
Definition: JDetector.hh:89
Router for direct addressing of module data in detector data structure.
Vec dir
track direction
Definition: Trk.hh:18
General purpose sorter of fit results.
*fatal Wrong number of arguments esac JCookie sh typeset Z DETECTOR typeset Z SOURCE_RUN typeset Z TARGET_RUN set_variable PARAMETERS_FILE $WORKDIR parameters
Definition: diff-Tuna.sh:38
Hit data structure.
General purpose class for parallel reading of objects from a single file or multiple files...
Empty structure for specification of parser element that is initialised (i.e. does not require input)...
Definition: JParser.hh:66
Long64_t counter_type
Type definition for counter.
double getTime(const Hit &hit)
Get true time of hit.
string outputFile
double E
Energy [GeV] (either MC truth or reconstructed)
Definition: Trk.hh:20
bool qualitySorter(const JFit &first, const JFit &second)
Comparison of fit results.
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
bool slewing
time slewing of analogue signal
The Vec class is a straightforward 3-d vector, which also works in pyroot.
Definition: Vec.hh:12
Detector file.
Definition: JHead.hh:196
JDirection3D getDirection(const Vec &dir)
Get direction.
Acoustic event fit.
int mc_run_id
MC run identifier.
Definition: Evt.hh:26
Auxiliary class for recursive type list generation.
Definition: JTypeList.hh:351
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1961
int getID() const
Get identifier.
Definition: JObjectID.hh:50
Data structure for PMT geometry, calibration and status.
Definition: JPMT.hh:43
JPosition3D getPosition(const Vec &pos)
Get position.
std::vector< double > fitinf
place to store additional fit info, for jgandalf, see JFitParameters.hh
Definition: Trk.hh:31
Auxiliary class for map of PMT parameters.
void load(const std::string &file_name)
Load mechanical model parameters from file.
Definition: JMechanics.hh:142
int debug
debug level
Definition: JSirene.cc:63
counter_type advance(counter_type &counter, const counter_type value, const counter_type limit=std::numeric_limits< counter_type >::max())
Advance counter.
Monte Carlo run header.
Definition: JHead.hh:1113
The Head class reflects the header of Monte-Carlo event files, which consists of keys (also referred ...
Definition: Head.hh:67
#define FATAL(A)
Definition: JMessage.hh:67
int id
track identifier
Definition: Trk.hh:16
Dynamic detector calibration.
Definition: JDynamics.hh:61
Vec pos
postion of the track at time t [m]
Definition: Trk.hh:17
Address of PMT in detector data structure.
Definition: JPMTAddress.hh:32
void load(const std::string &file_name, JDetector &detector)
Load detector from input file.
Data structure for set of track fit results.
std::vector< Trk > trks
list of reconstructed tracks (can be several because of prefits,showers, etc).
Definition: Evt.hh:36
double getLivetime(const std::string &file_name)
Get data taking live time.
std::vector< int > rec_stages
list of identifyers of succesfull fitting stages resulting in this track
Definition: Trk.hh:26
General purpose class for object reading from a list of file names.
std::vector< double > error_matrix
(NxN) error covariance matrix for fit parameters (stored as linear vector)
Definition: Trk.hh:33
double lik
likelihood or lambda value (for aafit, lambda)
Definition: Trk.hh:23
std::vector< Hit > mc_hits
MC: list of MC truth hits.
Definition: Evt.hh:45
static const int JPP_RECONSTRUCTION_TYPE
KM3NeT Data Definitions v2.1.0-7-g97845ea https://git.km3net.de/common/km3net-dataformat.
then usage $script< input_file >< detector_file > fi set_variable OUTPUT_DIR set_variable SELECTOR JDAQTimesliceL1 set_variable DEBUG case set_variable DEBUG
void copy(const Head &from, JHead &to)
Copy header from from to to.
Definition: JHead.cc:139
int id
offline event identifier
Definition: Evt.hh:21
Orientation of module.
const JLimit & getLimit() const
Get limit.
Definition: JLimit.hh:73
Data structure for PMT parameters.
std::vector< Hit > hits
list of hits
Definition: Evt.hh:35
do set_variable DETECTOR_TXT $WORKDIR detector
JAANET::DAQ DAQ
Definition: JHead.hh:1415
then fatal Wrong number of arguments fi set_variable DETECTOR $argv[1] set_variable INPUT_FILE $argv[2] eval JPrintDetector a $DETECTOR O IDENTIFIER eval JPrintDetector a $DETECTOR O SUMMARY source JAcoustics sh $DETECTOR_ID CHECK_EXIT_CODE typeset A TRIPODS get_tripods $WORKDIR tripod txt TRIPODS for EMITTER in
Definition: JCanberra.sh:41
The Trk class represents a Monte Carlo (MC) particle as well as a reconstructed track/shower.
Definition: Trk.hh:14
JTriggerCounter_t getCounter() const
Get trigger counter.
std::vector< Trk > mc_trks
MC: list of MC truth tracks.
Definition: Evt.hh:46
The Evt class respresent a Monte Carlo (MC) event as well as an offline event.
Definition: Evt.hh:19
then usage $script< detector file >< tripodfile >< stage > input file nInput files correspond to the output of JAcousticsEventBuilder[.sh] nFirst stage eval JPrintDetector a $DETECTOR O IDENTIFIER eval JPrintDetector a $DETECTOR O SUMMARY eval JPrintDetector a $DETECTOR O CAN source JAcoustics sh $DETECTOR_ID typeset A CONFIGURATION for key in Tmax_s
int rec_type
identifyer for the overall fitting algorithm/chain/strategy
Definition: Trk.hh:25