Fork me on GitHub

source: git/classes/DelphesSTDHEPReader.cc@ 84edab9

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 84edab9 was b55dc12, checked in by Pavel Demin <pavel.demin@…>, 10 years ago

simplify DelphesSTDHEPReader::SetInputFile

  • Property mode set to 100644
File size: 12.0 KB
Line 
1/*
2 * Delphes: a framework for fast simulation of a generic collider experiment
3 * Copyright (C) 2012-2014 Universite catholique de Louvain (UCL), Belgium
4 *
5 * This program is free software: you can redistribute it and/or modify
6 * it under the terms of the GNU General Public License as published by
7 * the Free Software Foundation, either version 3 of the License, or
8 * (at your option) any later version.
9 *
10 * This program is distributed in the hope that it will be useful,
11 * but WITHOUT ANY WARRANTY; without even the implied warranty of
12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13 * GNU General Public License for more details.
14 *
15 * You should have received a copy of the GNU General Public License
16 * along with this program. If not, see <http://www.gnu.org/licenses/>.
17 */
18
19
20/** \class DelphesSTDHEPReader
21 *
22 * Reads STDHEP file
23 *
24 * \author P. Demin - UCL, Louvain-la-Neuve
25 *
26 */
27
28#include "classes/DelphesSTDHEPReader.h"
29
30#include <stdexcept>
31#include <iostream>
32#include <sstream>
33
34#include <stdio.h>
35#include <errno.h>
36#include <rpc/types.h>
37#include <rpc/xdr.h>
38
39#include "TObjArray.h"
40#include "TStopwatch.h"
41#include "TDatabasePDG.h"
42#include "TParticlePDG.h"
43#include "TLorentzVector.h"
44
45#include "classes/DelphesClasses.h"
46#include "classes/DelphesFactory.h"
47
48#include "ExRootAnalysis/ExRootTreeBranch.h"
49
50using namespace std;
51
52static const int kBufferSize = 1000000;
53
54//---------------------------------------------------------------------------
55
56DelphesSTDHEPReader::DelphesSTDHEPReader() :
57 fInputFile(0), fInputXDR(0), fBuffer(0), fPDG(0), fBlockType(-1)
58{
59 fInputXDR = new XDR;
60 fBuffer = new char[kBufferSize*96 + 24];
61
62 fPDG = TDatabasePDG::Instance();
63}
64
65//---------------------------------------------------------------------------
66
67DelphesSTDHEPReader::~DelphesSTDHEPReader()
68{
69 if(fBuffer) delete fBuffer;
70 if(fInputXDR) delete fInputXDR;
71}
72
73//---------------------------------------------------------------------------
74
75void DelphesSTDHEPReader::SetInputFile(FILE *inputFile)
76{
77 fInputFile = inputFile;
78 xdrstdio_create(fInputXDR, inputFile, XDR_DECODE);
79}
80
81//---------------------------------------------------------------------------
82
83void DelphesSTDHEPReader::Clear()
84{
85 fBlockType = -1;
86}
87
88//---------------------------------------------------------------------------
89
90bool DelphesSTDHEPReader::EventReady()
91{
92 return (fBlockType == MCFIO_STDHEP) || (fBlockType == MCFIO_STDHEP4);
93}
94
95//---------------------------------------------------------------------------
96
97bool DelphesSTDHEPReader::ReadBlock(DelphesFactory *factory,
98 TObjArray *allParticleOutputArray,
99 TObjArray *stableParticleOutputArray,
100 TObjArray *partonOutputArray)
101{
102 if(feof(fInputFile)) return kFALSE;
103
104 xdr_int(fInputXDR, &fBlockType);
105
106 SkipBytes(4);
107
108 if(fBlockType == FILEHEADER)
109 {
110 ReadFileHeader();
111 }
112 else if(fBlockType == EVENTTABLE)
113 {
114 ReadEventTable();
115 }
116 else if(fBlockType == EVENTHEADER)
117 {
118 ReadEventHeader();
119 }
120 else if(fBlockType == MCFIO_STDHEPBEG ||
121 fBlockType == MCFIO_STDHEPEND)
122 {
123 ReadSTDCM1();
124 }
125 else if(fBlockType == MCFIO_STDHEP)
126 {
127 ReadSTDHEP();
128 AnalyzeParticles(factory, allParticleOutputArray,
129 stableParticleOutputArray, partonOutputArray);
130 }
131 else if(fBlockType == MCFIO_STDHEP4)
132 {
133 ReadSTDHEP();
134 AnalyzeParticles(factory, allParticleOutputArray,
135 stableParticleOutputArray, partonOutputArray);
136 ReadSTDHEP4();
137 }
138 else
139 {
140 throw runtime_error("Unsupported block type.");
141 }
142
143 return kTRUE;
144}
145
146//---------------------------------------------------------------------------
147
148void DelphesSTDHEPReader::SkipBytes(u_int size)
149{
150 int rc;
151 u_int rndup;
152
153 rndup = size % 4;
154 if(rndup > 0)
155 {
156 rndup = 4 - rndup;
157 }
158
159 rc = fseek(fInputFile, size + rndup, SEEK_CUR);
160
161 if(rc != 0 && errno == ESPIPE)
162 {
163 xdr_opaque(fInputXDR, fBuffer, size);
164 }
165}
166
167//---------------------------------------------------------------------------
168
169void DelphesSTDHEPReader::SkipArray(u_int elsize)
170{
171 u_int size;
172 xdr_u_int(fInputXDR, &size);
173 SkipBytes(size*elsize);
174}
175
176//---------------------------------------------------------------------------
177
178void DelphesSTDHEPReader::ReadFileHeader()
179{
180 u_int i;
181 enum STDHEPVersion {UNKNOWN, V1, V2, V21} version;
182
183 // version
184 xdr_string(fInputXDR, &fBuffer, 100);
185 if(fBuffer[0] == '\0' || fBuffer[1] == '\0') version = UNKNOWN;
186 else if(fBuffer[0] == '1') version = V1;
187 else if(strncmp(fBuffer, "2.01", 4) == 0) version = V21;
188 else if(fBuffer[0] == '2') version = V2;
189 else version = UNKNOWN;
190
191 if(version == UNKNOWN)
192 {
193 throw runtime_error("Unknown file format version.");
194 }
195
196 SkipArray(1);
197 SkipArray(1);
198 SkipArray(1);
199
200 if(version == V21)
201 {
202 SkipArray(1);
203 }
204
205 // Expected number of events
206 SkipBytes(4);
207
208 // Number of events
209 xdr_u_int(fInputXDR, &fEntries);
210
211 SkipBytes(8);
212
213 // Number of blocks
214 u_int nBlocks = 0;
215 xdr_u_int(fInputXDR, &nBlocks);
216
217 // Number of NTuples
218 u_int nNTuples = 0;
219 if(version != V1)
220 {
221 xdr_u_int(fInputXDR, &nNTuples);
222 }
223
224 if(nNTuples != 0)
225 {
226 throw runtime_error("Files containing n-tuples are not supported.");
227 }
228
229 // Processing blocks extraction
230 if(nBlocks != 0)
231 {
232 SkipArray(4);
233
234 for(i = 0; i < nBlocks; i++)
235 {
236 SkipArray(1);
237 }
238 }
239}
240
241//---------------------------------------------------------------------------
242
243void DelphesSTDHEPReader::ReadEventTable()
244{
245 // version
246 xdr_string(fInputXDR, &fBuffer, 100);
247 if(strncmp(fBuffer, "1.00", 4) == 0)
248 {
249 SkipBytes(8);
250
251 SkipArray(4);
252 SkipArray(4);
253 SkipArray(4);
254 SkipArray(4);
255 SkipArray(4);
256 }
257 else if(strncmp(fBuffer, "2.00", 4) == 0)
258 {
259 SkipBytes(12);
260
261 SkipArray(4);
262 SkipArray(4);
263 SkipArray(4);
264 SkipArray(4);
265 SkipArray(8);
266 }
267}
268
269//---------------------------------------------------------------------------
270
271void DelphesSTDHEPReader::ReadEventHeader()
272{
273 bool skipNTuples = false;
274 u_int skipSize = 4;
275
276 // version
277 xdr_string(fInputXDR, &fBuffer, 100);
278 if(strncmp(fBuffer, "2.00", 4) == 0)
279 {
280 skipNTuples = true;
281 }
282 else if(strncmp(fBuffer, "3.00", 4) == 0)
283 {
284 skipNTuples = true;
285 skipSize = 8;
286 }
287
288 SkipBytes(20);
289
290 u_int dimBlocks = 0;
291 xdr_u_int(fInputXDR, &dimBlocks);
292
293 u_int dimNTuples = 0;
294 if(skipNTuples)
295 {
296 SkipBytes(4);
297 xdr_u_int(fInputXDR, &dimNTuples);
298 }
299
300 // Processing blocks extraction
301 if(dimBlocks > 0)
302 {
303 SkipArray(4);
304 SkipArray(skipSize);
305 }
306
307 // Processing blocks extraction
308 if(skipNTuples && dimNTuples > 0)
309 {
310 SkipArray(4);
311 SkipArray(skipSize);
312 }
313}
314
315//---------------------------------------------------------------------------
316
317void DelphesSTDHEPReader::ReadSTDCM1()
318{
319 // version
320 xdr_string(fInputXDR, &fBuffer, 100);
321
322 // skip 5*4 + 2*8 = 36 bytes
323 SkipBytes(36);
324
325 if((strncmp(fBuffer, "1.", 2) == 0) || (strncmp(fBuffer, "2.", 2) == 0) ||
326 (strncmp(fBuffer, "3.", 2) == 0) || (strncmp(fBuffer, "4.", 2) == 0) ||
327 (strncmp(fBuffer, "5.00", 4) == 0))
328 {
329 return;
330 }
331
332 SkipArray(1);
333 SkipArray(1);
334
335 if(strncmp(fBuffer, "5.01", 4) == 0)
336 {
337 return;
338 }
339
340 SkipBytes(4);
341}
342
343//---------------------------------------------------------------------------
344
345void DelphesSTDHEPReader::ReadSTDHEP()
346{
347 u_int idhepSize, isthepSize, jmohepSize, jdahepSize, phepSize, vhepSize;
348
349 // version
350 xdr_string(fInputXDR, &fBuffer, 100);
351
352 // Extracting the event number
353 xdr_int(fInputXDR, &fEventNumber);
354
355 // Extracting the number of particles
356 xdr_int(fInputXDR, &fEventSize);
357
358 if(fEventSize >= kBufferSize)
359 {
360 throw runtime_error("too many particles in event");
361 }
362
363 // 4*n + 4*n + 8*n + 8*n + 40*n + 32*n +
364 // 4 + 4 + 4 + 4 + 4 + 4 = 96*n + 24
365
366 xdr_opaque(fInputXDR, fBuffer, 96*fEventSize + 24);
367
368 idhepSize = ntohl(*(u_int*)(fBuffer));
369 isthepSize = ntohl(*(u_int*)(fBuffer + 4*1 + 4*1*fEventSize));
370 jmohepSize = ntohl(*(u_int*)(fBuffer + 4*2 + 4*2*fEventSize));
371 jdahepSize = ntohl(*(u_int*)(fBuffer + 4*3 + 4*4*fEventSize));
372 phepSize = ntohl(*(u_int*)(fBuffer + 4*4 + 4*6*fEventSize));
373 vhepSize = ntohl(*(u_int*)(fBuffer + 4*5 + 4*16*fEventSize));
374
375 if(fEventSize < 0 ||
376 fEventSize != (int)idhepSize || fEventSize != (int)isthepSize ||
377 (2*fEventSize) != (int)jmohepSize || (2*fEventSize) != (int)jdahepSize ||
378 (5*fEventSize) != (int)phepSize || (4*fEventSize) != (int)vhepSize)
379 {
380 throw runtime_error("Inconsistent size of arrays. File is probably corrupted.");
381 }
382
383 fWeight = 1.0;
384 fAlphaQED = 0.0;
385 fAlphaQCD = 0.0;
386 fScaleSize = 0;
387 memset(fScale, 0, 10*sizeof(double));
388}
389
390//---------------------------------------------------------------------------
391
392void DelphesSTDHEPReader::ReadSTDHEP4()
393{
394 u_int number;
395
396 // Extracting the event weight
397 xdr_double(fInputXDR, &fWeight);
398
399 // Extracting alpha QED
400 xdr_double(fInputXDR, &fAlphaQED);
401
402 // Extracting alpha QCD
403 xdr_double(fInputXDR, &fAlphaQCD);
404
405 // Extracting the event scale
406 xdr_u_int(fInputXDR, &fScaleSize);
407 for(number = 0; number < fScaleSize; ++number)
408 {
409 xdr_double(fInputXDR, &fScale[number]);
410 }
411
412 SkipArray(8);
413 SkipArray(4);
414
415 SkipBytes(4);
416}
417
418//---------------------------------------------------------------------------
419
420void DelphesSTDHEPReader::AnalyzeEvent(ExRootTreeBranch *branch, long long eventNumber,
421 TStopwatch *readStopWatch, TStopwatch *procStopWatch)
422{
423 LHEFEvent *element;
424
425 element = static_cast<LHEFEvent *>(branch->NewEntry());
426
427 element->Number = fEventNumber;
428
429 element->ProcessID = 0;
430
431 element->Weight = fWeight;
432 element->ScalePDF = fScale[0];
433 element->AlphaQED = fAlphaQED;
434 element->AlphaQCD = fAlphaQCD;
435
436 element->ReadTime = readStopWatch->RealTime();
437 element->ProcTime = procStopWatch->RealTime();
438}
439
440//---------------------------------------------------------------------------
441
442void DelphesSTDHEPReader::AnalyzeParticles(DelphesFactory *factory,
443 TObjArray *allParticleOutputArray,
444 TObjArray *stableParticleOutputArray,
445 TObjArray *partonOutputArray)
446{
447 Candidate *candidate;
448 TParticlePDG *pdgParticle;
449 int pdgCode;
450
451 int number;
452 int pid, status, m1, m2, d1, d2;
453 double px, py, pz, e, mass;
454 double x, y, z, t;
455
456 XDR bufferXDR[6];
457 xdrmem_create(&bufferXDR[0], fBuffer + 4*1, 4*fEventSize, XDR_DECODE);
458 xdrmem_create(&bufferXDR[1], fBuffer + 4*2 + 4*1*fEventSize, 4*fEventSize, XDR_DECODE);
459 xdrmem_create(&bufferXDR[2], fBuffer + 4*3 + 4*2*fEventSize, 8*fEventSize, XDR_DECODE);
460 xdrmem_create(&bufferXDR[3], fBuffer + 4*4 + 4*4*fEventSize, 8*fEventSize, XDR_DECODE);
461 xdrmem_create(&bufferXDR[4], fBuffer + 4*5 + 4*6*fEventSize, 40*fEventSize, XDR_DECODE);
462 xdrmem_create(&bufferXDR[5], fBuffer + 4*6 + 4*16*fEventSize, 32*fEventSize, XDR_DECODE);
463
464 for(number = 0; number < fEventSize; ++number)
465 {
466 xdr_int(&bufferXDR[0], &status);
467 xdr_int(&bufferXDR[1], &pid);
468 xdr_int(&bufferXDR[2], &m1);
469 xdr_int(&bufferXDR[2], &m2);
470 xdr_int(&bufferXDR[3], &d1);
471 xdr_int(&bufferXDR[3], &d2);
472
473 xdr_double(&bufferXDR[4], &px);
474 xdr_double(&bufferXDR[4], &py);
475 xdr_double(&bufferXDR[4], &pz);
476 xdr_double(&bufferXDR[4], &e);
477 xdr_double(&bufferXDR[4], &mass);
478
479 xdr_double(&bufferXDR[5], &x);
480 xdr_double(&bufferXDR[5], &y);
481 xdr_double(&bufferXDR[5], &z);
482 xdr_double(&bufferXDR[5], &t);
483
484 candidate = factory->NewCandidate();
485
486 candidate->PID = pid;
487 pdgCode = TMath::Abs(candidate->PID);
488
489 candidate->Status = status;
490
491 candidate->M1 = m1 - 1;
492 candidate->M2 = m2 - 1;
493
494 candidate->D1 = d1 - 1;
495 candidate->D2 = d2 - 1;
496
497 pdgParticle = fPDG->GetParticle(pid);
498 candidate->Charge = pdgParticle ? int(pdgParticle->Charge()/3.0) : -999;
499 candidate->Mass = mass;
500
501 candidate->Momentum.SetPxPyPzE(px, py, pz, e);
502
503 candidate->Position.SetXYZT(x, y, z, t);
504
505 allParticleOutputArray->Add(candidate);
506
507 if(!pdgParticle) continue;
508
509 if(status == 1 && pdgParticle->Stable())
510 {
511 stableParticleOutputArray->Add(candidate);
512 }
513 else if(pdgCode <= 5 || pdgCode == 21 || pdgCode == 15)
514 {
515 partonOutputArray->Add(candidate);
516 }
517 }
518}
519
520//---------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.