Jpp 19.3.0-rc.1
the software that should make you happy
Loading...
Searching...
No Matches
JDST.cc File Reference

Example program to select events from DST files. More...

#include <cstring>
#include <string>
#include <iostream>
#include <iomanip>
#include <vector>
#include <set>
#include "TError.h"
#include "TROOT.h"
#include "TFile.h"
#include "TKey.h"
#include "TTree.h"
#include "TChain.h"
#include "km3net-dataformat/offline/Evt.hh"
#include "km3net-dataformat/definitions/root.hh"
#include "JSystem/JGlob.hh"
#include "JSupport/JLimit.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

Example program to select events from DST files.

Author
mdejong

Definition in file JDST.cc.

Function Documentation

◆ main()

int main ( int argc,
char ** argv )

Definition at line 120 of file JDST.cc.

121{
122 using namespace std;
123 using namespace JPP;
124
125 vector<string> inputFile;
126 JLimit_t numberOfEvents;
127 string selectFile;
128 string outputFile;
129 int debug;
130
131 try {
132
133 JParser<> zap("Example program to select events from DST files.");
134
135 zap['f'] = make_field(inputFile);
136 zap['n'] = make_field(numberOfEvents) = JLimit::max();
137 zap['<'] = make_field(selectFile,
138 "<run> <frame index> <trigger counter>"\
139 " (one header line)") = "";
140 zap['o'] = make_field(outputFile);
141 zap['d'] = make_field(debug) = 1;
142
143 zap(argc, argv);
144 }
145 catch(const exception &error) {
146 FATAL(error.what() << endl);
147 }
148
149 gErrorIgnoreLevel = kFatal;
150
151 inputFile = getFilenames(inputFile);
152
153 set<JTag> selection;
154
155 if (selectFile != "") {
156
157 ifstream in(selectFile.c_str());
158
159 in.ignore(numeric_limits<streamsize>::max(), '\n');
160
161 for (JTag tag; in >> tag; ) {
162 selection.insert(tag);
163 }
164
165 in.close();
166 }
167
168 TFile* out = TFile::Open(outputFile.c_str(), "CREATE");
169
170 if (out == NULL || ! out->IsOpen()) {
171 FATAL("File " << outputFile << " not opened." << endl);
172 }
173
174 vector<TChain*> input;
175 vector<TTree*> output;
176
177 for (vector<string>::const_iterator file_name = inputFile.begin(); file_name != inputFile.end(); ++file_name) {
178
179 TFile* in = TFile::Open(file_name->c_str(), "EXISTS");
180
181 if (in != NULL && in->IsOpen()) {
182
183 TIter iter(in->GetListOfKeys(), kIterBackward);
184
185 for (TKey* key; (key = (TKey*) iter.Next()) != NULL; ) {
186
187 TKey* p = dynamic_cast<TKey*>(in->GetListOfKeys()->Before(key));
188
189 if (p == NULL || strcmp(key->GetName(), p->GetName()) != 0) { // select last key
190
191 TTree* t1 = dynamic_cast<TTree*>(key->ReadObj());
192
193 if (t1 != NULL) {
194
195 TChain* q = NULL;
196
197 for (vector<TChain*>::iterator i = input.begin(); i != input.end(); ++i) {
198
199 if (strcmp((*i)->GetName(), t1->GetName()) == 0) {
200
201 q = *i;
202
203 break;
204 }
205 }
206
207 if (q == NULL) {
208 input.push_back(q = new TChain(t1->GetName()));
209 }
210
211 q->Add(file_name->c_str());
212 }
213 }
214 }
215 }
216 }
217
218 // master
219
220 TChain* tm = NULL;
221
222 for (vector<TChain*>::iterator i = input.begin(); i != input.end(); ++i) {
223
224 if (strcmp((*i)->GetName(), TTREE_OFFLINE_EVENT) == 0) {
225
226 tm = *i;
227
228 break;
229 }
230 }
231
232 if (tm == NULL) {
233 FATAL("Missing TTree " << TTREE_OFFLINE_EVENT << endl);
234 }
235
236 Evt* evt = NULL;
237
238 tm->GetEntries(); // load TTree
239 tm->SetBranchAddress(TBRANCH_OFFLINE_EVENT, &evt); //
240
241 // select friend TTree's based on number of entries
242
243 for (vector<TChain*>::iterator i = input.begin(); i != input.end(); ) {
244 if (tm->GetEntries() == (*i)->GetEntries())
245 ++i;
246 else
247 i = input.erase(i);
248 }
249
250 out->cd();
251
252 for (vector<TChain*>::iterator i = input.begin(); i != input.end(); ++i) {
253 output.push_back((*i)->GetTree()->CloneTree(0));
254 }
255
256 for (Long64_t i = numberOfEvents.getLowerLimit(); i < tm->GetEntries() && i < numberOfEvents.getUpperLimit(); ++i) {
257
258 STATUS("event: " << setw(10) << i << '\r'); DEBUG(endl);
259
260 tm->GetEntry(i);
261
262 const JTag tag(evt->run_id, evt->frame_index, evt->trigger_counter);
263
264 DEBUG("tag " << tag << ' ' << (selection.empty() || selection.count(tag)) << endl);
265
266 if (selection.empty() || selection.count(tag)) {
267
268 for (size_t m = 0; m != input.size(); ++m) {
269 input [m]->GetEntry(i);
270 output[m]->Fill();
271 }
272 }
273 }
274 STATUS(endl);
275
276 out->Write();
277 out->Close();
278}
string outputFile
#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
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition JParser.hh:2142
Utility class to parse command line options.
Definition JParser.hh:1698
T getLowerLimit() const
Get lower limit.
Definition JRange.hh:202
T getUpperLimit() const
Get upper limit.
Definition JRange.hh:213
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
static const char *const TTREE_OFFLINE_EVENT
ROOT TTree name.
Definition root.hh:19
static const char *const TBRANCH_OFFLINE_EVENT
ROOT TBranch name.
Definition root.hh:28
The Evt class respresent a Monte Carlo (MC) event as well as an offline event.
Definition Evt.hh:21
int frame_index
from the raw data
Definition Evt.hh:29
int run_id
DAQ run identifier.
Definition Evt.hh:26
ULong64_t trigger_counter
trigger counter
Definition Evt.hh:31
Auxiliary class for defining the range of iterations of objects.
Definition JLimit.hh:45