Jpp
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Functions
JMergeFit.cc File Reference

Program to merge different files with JFIT::JEvt data. More...

#include <string>
#include <iostream>
#include <iomanip>
#include <vector>
#include <algorithm>
#include <iterator>
#include "km3net-dataformat/offline/Head.hh"
#include "km3net-dataformat/offline/Evt.hh"
#include "JDAQ/JDAQEventIO.hh"
#include "JDAQ/JDAQTimesliceIO.hh"
#include "JDAQ/JDAQSummarysliceIO.hh"
#include "JDAQ/JDAQEvaluator.hh"
#include "JTrigger/JTriggerParameters.hh"
#include "JSupport/JParallelFileScanner.hh"
#include "JSupport/JTreeScanner.hh"
#include "JSupport/JFileRecorder.hh"
#include "JSupport/JSupport.hh"
#include "JSupport/JMeta.hh"
#include "JReconstruction/JEvt.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

Program to merge different files with JFIT::JEvt data.

Author
mdejong

Definition in file JMergeFit.cc.

Function Documentation

int main ( int  argc,
char **  argv 
)

Definition at line 87 of file JMergeFit.cc.

88 {
89  using namespace std;
90  using namespace JPP;
91  using namespace KM3NETDAQ;
92 
93  typedef JParallelFileScanner< JTypeList<JDAQEvent, JEvt> > JParallelFileScanner_t;
94  typedef JParallelFileScanner_t::multi_pointer_type multi_pointer_type;
96 
97  vector<std::string> inputFile;
99  int debug;
100 
101  try {
102 
103  JParser<> zap("Program to merge different files with JFIT::JEvt data.");
104 
105  zap['f'] = make_field(inputFile, "list of JEvt compatible files");
106  zap['o'] = make_field(outputFile, "single file with merged JEvt data");
107  zap['d'] = make_field(debug) = 1;
108 
109  zap(argc, argv);
110  }
111  catch(const exception& error) {
112  FATAL(error.what() << endl);
113  }
114 
115  cout.tie(&cerr);
116 
117  if (inputFile.empty()) {
118  FATAL("No input files." << endl);
119  }
120 
121  outputFile.open();
122 
123  outputFile.put(JMeta(argc, argv));
124 
125 
126  JParallelFileScanner_t in(inputFile[0]);
127 
128  vector<JTreeRouter_t> buffer(1); // skip one entry
129 
130  for (size_t i = 1; i != inputFile.size(); ++i) {
131  buffer.push_back(JTreeRouter_t(inputFile[i]));
132  }
133 
134  while (in.hasNext()) {
135 
136  STATUS("event: " << setw(10) << in.getCounter() << '\r'); DEBUG(endl);
137 
138  multi_pointer_type ps = in.next();
139 
140  const JDAQEvent* tev = ps;
141  const JEvt* evt = ps;
142 
143  JEvt out(*evt);
144 
145  for (size_t i = 1; i != inputFile.size(); ++i) {
146 
147  const JEvt* p = buffer[i].get(*tev);
148 
149  if (p != NULL) {
150 
151  copy(p->begin(), p->end(), back_inserter(out));
152 
153  } else {
154 
155  FATAL("Inconsistent data at " << in.getFilename() << ":" << in.getCounter() << endl);
156  }
157  }
158 
159  outputFile.put(out);
160  }
161  STATUS(endl);
162 
164 
165  io >> outputFile;
166 
167  outputFile.close();
168 }
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
#define STATUS(A)
Definition: JMessage.hh:63
General purpose class for parallel reading of objects from a single file or multiple files...
string outputFile
Type list.
Definition: JTypeList.hh:22
Acoustic event fit.
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1961
int debug
debug level
Definition: JSirene.cc:63
#define FATAL(A)
Definition: JMessage.hh:67
General purpose class for object reading from a list of file names.
void copy(const Head &from, JHead &to)
Copy header from from to to.
Definition: JHead.cc:139
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 typeset A TRIPODS get_tripods $WORKDIR tripod txt TRIPODS for EMITTER in
Definition: JCanberra.sh:36
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62