Fork me on GitHub

source: git/classes/DelphesHepMC3Reader.cc@ 26a21bd

Last change on this file since 26a21bd was a26a130, checked in by Pavel Demin <pavel.demin@…>, 3 years ago

move arrays to FinalizeParticles in DelphesHepMC3Reader

  • Property mode set to 100644
File size: 13.5 KB
RevLine 
[b443089]1/*
2 * Delphes: a framework for fast simulation of a generic collider experiment
[b4786d3]3 * Copyright (C) 2012-2021 Universite catholique de Louvain (UCL), Belgium
[1fa50c2]4 *
[b443089]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.
[1fa50c2]9 *
[b443089]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.
[1fa50c2]14 *
[b443089]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
[b4786d3]19/** \class DelphesHepMC3Reader
[d7d2da3]20 *
21 * Reads HepMC file
22 *
23 * \author P. Demin - UCL, Louvain-la-Neuve
24 *
25 */
26
[b4786d3]27#include "classes/DelphesHepMC3Reader.h"
[d7d2da3]28
29#include <iostream>
30#include <sstream>
[341014c]31#include <stdexcept>
[d7d2da3]32
33#include <map>
[59abd43]34#include <vector>
[d7d2da3]35
36#include <stdio.h>
37
38#include "TDatabasePDG.h"
39#include "TLorentzVector.h"
[341014c]40#include "TObjArray.h"
41#include "TParticlePDG.h"
42#include "TStopwatch.h"
[d7d2da3]43
44#include "classes/DelphesClasses.h"
45#include "classes/DelphesFactory.h"
46#include "classes/DelphesStream.h"
47
48#include "ExRootAnalysis/ExRootTreeBranch.h"
49
50using namespace std;
51
[c6667c0]52static const int kBufferSize = 16384;
[d7d2da3]53
54//---------------------------------------------------------------------------
55
[b4786d3]56DelphesHepMC3Reader::DelphesHepMC3Reader() :
[d7d2da3]57 fInputFile(0), fBuffer(0), fPDG(0),
[4006893]58 fVertexCounter(-2), fParticleCounter(-1)
[d7d2da3]59{
60 fBuffer = new char[kBufferSize];
61
62 fPDG = TDatabasePDG::Instance();
63}
64
65//---------------------------------------------------------------------------
66
[b4786d3]67DelphesHepMC3Reader::~DelphesHepMC3Reader()
[d7d2da3]68{
69 if(fBuffer) delete[] fBuffer;
70}
71
72//---------------------------------------------------------------------------
73
[b4786d3]74void DelphesHepMC3Reader::SetInputFile(FILE *inputFile)
[d7d2da3]75{
76 fInputFile = inputFile;
77}
78
79//---------------------------------------------------------------------------
80
[b4786d3]81void DelphesHepMC3Reader::Clear()
[d7d2da3]82{
[4006893]83 fWeights.clear();
[cf22deb]84 fMomentumCoefficient = 1.0;
85 fPositionCoefficient = 1.0;
[4006893]86 fVertexCounter = -2;
[b4786d3]87 fParticleCounter = -1;
[4006893]88 fVertices.clear();
89 fParticles.clear();
90 fInVertexMap.clear();
91 fOutVertexMap.clear();
92 fMotherMap.clear();
[d7d2da3]93 fDaughterMap.clear();
94}
95
96//---------------------------------------------------------------------------
97
[b4786d3]98bool DelphesHepMC3Reader::EventReady()
[d7d2da3]99{
[4006893]100 return (fVertexCounter == -1) && (fParticleCounter == 0);
[d7d2da3]101}
102
103//---------------------------------------------------------------------------
104
[b4786d3]105bool DelphesHepMC3Reader::ReadBlock(DelphesFactory *factory,
[d7d2da3]106 TObjArray *allParticleOutputArray,
107 TObjArray *stableParticleOutputArray,
108 TObjArray *partonOutputArray)
109{
[77e9ae1]110 map<int, pair<int, int> >::iterator itDaughterMap;
[cf22deb]111 char key, momentumUnit[4], positionUnit[3];
[b4786d3]112 int rc, code;
[59abd43]113 double weight;
[d7d2da3]114
115 if(!fgets(fBuffer, kBufferSize, fInputFile)) return kFALSE;
116
117 DelphesStream bufferStream(fBuffer + 1);
118
119 key = fBuffer[0];
120
121 if(key == 'E')
122 {
123 Clear();
124
125 rc = bufferStream.ReadInt(fEventNumber)
[59abd43]126 && bufferStream.ReadInt(fVertexCounter)
[b4786d3]127 && bufferStream.ReadInt(fParticleCounter);
[59abd43]128
129 if(!rc)
130 {
[341014c]131 cerr << "** ERROR: "
132 << "invalid event format" << endl;
[59abd43]133 return kFALSE;
134 }
[b4786d3]135 }
136 else if(key == 'U')
137 {
138 rc = sscanf(fBuffer + 1, "%3s %2s", momentumUnit, positionUnit);
[59abd43]139
[b4786d3]140 if(rc != 2)
[59abd43]141 {
[b4786d3]142 cerr << "** ERROR: "
143 << "invalid units format" << endl;
144 return kFALSE;
[59abd43]145 }
146
[b4786d3]147 if(strncmp(momentumUnit, "GEV", 3) == 0)
[59abd43]148 {
[b4786d3]149 fMomentumCoefficient = 1.0;
150 }
151 else if(strncmp(momentumUnit, "MEV", 3) == 0)
152 {
153 fMomentumCoefficient = 0.001;
[59abd43]154 }
155
[b4786d3]156 if(strncmp(positionUnit, "MM", 3) == 0)
157 {
158 fPositionCoefficient = 1.0;
159 }
160 else if(strncmp(positionUnit, "CM", 3) == 0)
161 {
162 fPositionCoefficient = 10.0;
163 }
164 }
165 else if(key == 'W')
166 {
167 while(bufferStream.ReadDbl(weight))
[59abd43]168 {
[4006893]169 fWeights.push_back(weight);
[59abd43]170 }
[b4786d3]171 }
172 else if(key == 'A' && bufferStream.FindStr("mpi"))
173 {
174 rc = bufferStream.ReadInt(fMPI);
[d7d2da3]175
176 if(!rc)
177 {
[341014c]178 cerr << "** ERROR: "
[b4786d3]179 << "invalid MPI format" << endl;
[d7d2da3]180 return kFALSE;
181 }
182 }
[b4786d3]183 else if(key == 'A' && bufferStream.FindStr("signal_process_id"))
[cf22deb]184 {
[b4786d3]185 rc = bufferStream.ReadInt(fProcessID);
[cf22deb]186
[b4786d3]187 if(!rc)
[cf22deb]188 {
[341014c]189 cerr << "** ERROR: "
[b4786d3]190 << "invalid process ID format" << endl;
[cf22deb]191 return kFALSE;
192 }
[b4786d3]193 }
194 else if(key == 'A' && bufferStream.FindStr("event_scale"))
195 {
196 rc = bufferStream.ReadDbl(fScale);
[cf22deb]197
[b4786d3]198 if(!rc)
[cf22deb]199 {
[b4786d3]200 cerr << "** ERROR: "
201 << "invalid event scale format" << endl;
202 return kFALSE;
[cf22deb]203 }
[b4786d3]204 }
205 else if(key == 'A' && bufferStream.FindStr("alphaQCD"))
206 {
207 rc = bufferStream.ReadDbl(fAlphaQCD);
[341014c]208
[b4786d3]209 if(!rc)
[cf22deb]210 {
[b4786d3]211 cerr << "** ERROR: "
212 << "invalid alphaQCD format" << endl;
213 return kFALSE;
[cf22deb]214 }
[b4786d3]215 }
216 else if(key == 'A' && bufferStream.FindStr("alphaQED"))
217 {
218 rc = bufferStream.ReadDbl(fAlphaQED);
219
220 if(!rc)
[cf22deb]221 {
[b4786d3]222 cerr << "** ERROR: "
223 << "invalid alphaQED format" << endl;
224 return kFALSE;
[cf22deb]225 }
226 }
[b4786d3]227 else if(key == 'A' && bufferStream.FindStr("GenCrossSection"))
[edeb0f0]228 {
229 rc = bufferStream.ReadDbl(fCrossSection)
230 && bufferStream.ReadDbl(fCrossSectionError);
[341014c]231
[b4786d3]232 if(!rc)
233 {
234 cerr << "** ERROR: "
235 << "invalid cross section format" << endl;
236 return kFALSE;
237 }
238 }
239 else if(key == 'A' && bufferStream.FindStr("GenPdfInfo"))
[d7d2da3]240 {
241 rc = bufferStream.ReadInt(fID1)
242 && bufferStream.ReadInt(fID2)
243 && bufferStream.ReadDbl(fX1)
244 && bufferStream.ReadDbl(fX2)
245 && bufferStream.ReadDbl(fScalePDF)
246 && bufferStream.ReadDbl(fPDF1)
247 && bufferStream.ReadDbl(fPDF2);
248
249 if(!rc)
250 {
[341014c]251 cerr << "** ERROR: "
252 << "invalid PDF format" << endl;
[d7d2da3]253 return kFALSE;
254 }
255 }
[b4786d3]256 else if(key == 'V')
[d7d2da3]257 {
[4006893]258 fParticles.clear();
259
[b4786d3]260 fX = 0.0;
261 fY = 0.0;
262 fZ = 0.0;
263 fT = 0.0;
264
265 rc = bufferStream.ReadInt(fVertexCode)
266 && bufferStream.ReadInt(fVertexStatus);
267
268 if(!rc)
269 {
270 cerr << "** ERROR: "
271 << "invalid vertex format" << endl;
272 return kFALSE;
273 }
274
275 rc = bufferStream.FindChr('[');
[d7d2da3]276
277 if(!rc)
278 {
[341014c]279 cerr << "** ERROR: "
280 << "invalid vertex format" << endl;
[d7d2da3]281 return kFALSE;
282 }
[b4786d3]283
284 while(bufferStream.ReadInt(code))
285 {
[4006893]286 fParticles.push_back(code);
287 bufferStream.FindChr(',');
[b4786d3]288 }
289
290 if(bufferStream.FindChr('@'))
291 {
292 rc = bufferStream.ReadDbl(fX)
293 && bufferStream.ReadDbl(fY)
294 && bufferStream.ReadDbl(fZ)
295 && bufferStream.ReadDbl(fT);
296
297 if(!rc)
298 {
299 cerr << "** ERROR: "
300 << "invalid vertex format" << endl;
301 return kFALSE;
302 }
303 }
[4006893]304
305 AnalyzeVertex(factory, fVertexCode);
[d7d2da3]306 }
[b4786d3]307 else if(key == 'P' && fParticleCounter > 0)
[d7d2da3]308 {
[b4786d3]309 --fParticleCounter;
310
[d7d2da3]311 rc = bufferStream.ReadInt(fParticleCode)
[b4786d3]312 && bufferStream.ReadInt(fOutVertexCode)
[d7d2da3]313 && bufferStream.ReadInt(fPID)
314 && bufferStream.ReadDbl(fPx)
315 && bufferStream.ReadDbl(fPy)
316 && bufferStream.ReadDbl(fPz)
317 && bufferStream.ReadDbl(fE)
318 && bufferStream.ReadDbl(fMass)
[b4786d3]319 && bufferStream.ReadInt(fParticleStatus);
[d7d2da3]320
321 if(!rc)
322 {
[341014c]323 cerr << "** ERROR: "
324 << "invalid particle format" << endl;
[d7d2da3]325 return kFALSE;
326 }
327
[a26a130]328 AnalyzeParticle(factory);
[d7d2da3]329 }
330
331 if(EventReady())
332 {
[a26a130]333 FinalizeParticles(allParticleOutputArray, stableParticleOutputArray, partonOutputArray);
[d7d2da3]334 }
335
336 return kTRUE;
337}
338
339//---------------------------------------------------------------------------
340
[b4786d3]341void DelphesHepMC3Reader::AnalyzeEvent(ExRootTreeBranch *branch, long long eventNumber,
[d7d2da3]342 TStopwatch *readStopWatch, TStopwatch *procStopWatch)
343{
344 HepMCEvent *element;
345
346 element = static_cast<HepMCEvent *>(branch->NewEntry());
347 element->Number = fEventNumber;
348
349 element->ProcessID = fProcessID;
350 element->MPI = fMPI;
[4006893]351 element->Weight = fWeights.size() > 0 ? fWeights[0] : 1.0;
[edeb0f0]352 element->CrossSection = fCrossSection;
353 element->CrossSectionError = fCrossSectionError;
[d7d2da3]354 element->Scale = fScale;
355 element->AlphaQED = fAlphaQED;
356 element->AlphaQCD = fAlphaQCD;
357
358 element->ID1 = fID1;
359 element->ID2 = fID2;
360 element->X1 = fX1;
361 element->X2 = fX2;
362 element->ScalePDF = fScalePDF;
363 element->PDF1 = fPDF1;
364 element->PDF2 = fPDF2;
365
366 element->ReadTime = readStopWatch->RealTime();
367 element->ProcTime = procStopWatch->RealTime();
368}
369
370//---------------------------------------------------------------------------
371
[b4786d3]372void DelphesHepMC3Reader::AnalyzeWeight(ExRootTreeBranch *branch)
[d203a37]373{
374 Weight *element;
[a26a130]375 vector<double>::const_iterator itWeights;
[d203a37]376
[a26a130]377 for(itWeights = fWeights.begin(); itWeights != fWeights.end(); ++itWeights)
[d203a37]378 {
379 element = static_cast<Weight *>(branch->NewEntry());
380
[a26a130]381 element->Weight = *itWeights;
[d203a37]382 }
383}
384
385//---------------------------------------------------------------------------
386
[4006893]387void DelphesHepMC3Reader::AnalyzeVertex(DelphesFactory *factory, int code, Candidate *candidate)
388{
389 int index;
390 TLorentzVector *position;
391 TObjArray *array;
392 vector<int>::iterator itParticle;
393 map<int, int>::iterator itVertexMap;
394
395 itVertexMap = fOutVertexMap.find(code);
396 if(itVertexMap == fOutVertexMap.end())
397 {
398 --fVertexCounter;
399
400 index = fVertices.size();
401 fOutVertexMap[code] = index;
402 if(candidate && code > 0) fInVertexMap[code] = index;
403
404 position = factory->New<TLorentzVector>();
405 array = factory->NewArray();
406 position->SetXYZT(0.0, 0.0, 0.0, 0.0);
407 fVertices.push_back(make_pair(position, array));
408 }
409 else
410 {
411 index = itVertexMap->second;
412 position = fVertices[index].first;
413 array = fVertices[index].second;
414 }
415
416 if(candidate)
417 {
418 array->Add(candidate);
419 }
420 else
421 {
422 position->SetXYZT(fX, fY, fZ, fT);
423 for(itParticle = fParticles.begin(); itParticle != fParticles.end(); ++itParticle)
424 {
425 fInVertexMap[*itParticle] = index;
426 }
427 }
428}
429
430//---------------------------------------------------------------------------
431
[a26a130]432void DelphesHepMC3Reader::AnalyzeParticle(DelphesFactory *factory)
[d7d2da3]433{
434 Candidate *candidate;
435
436 candidate = factory->NewCandidate();
437
438 candidate->PID = fPID;
439
[b4786d3]440 candidate->Status = fParticleStatus;
[d7d2da3]441
[80d4a34]442 candidate->Mass = fMass;
[d7d2da3]443
444 candidate->Momentum.SetPxPyPzE(fPx, fPy, fPz, fE);
445
[4006893]446 candidate->D1 = fParticleCode;
[b4786d3]447
[4006893]448 AnalyzeVertex(factory, fOutVertexCode, candidate);
[d7d2da3]449}
450
451//---------------------------------------------------------------------------
452
[a26a130]453void DelphesHepMC3Reader::FinalizeParticles(TObjArray *allParticleOutputArray,
454 TObjArray *stableParticleOutputArray,
455 TObjArray *partonOutputArray)
[d7d2da3]456{
[4006893]457 TLorentzVector *position;
458 TObjArray *array;
[d7d2da3]459 Candidate *candidate;
[a26a130]460 TParticlePDG *pdgParticle;
461 int pdgCode;
[b4786d3]462 map<int, int >::iterator itVertexMap;
[4006893]463 map<int, pair<int, int> >::iterator itMotherMap;
[77e9ae1]464 map<int, pair<int, int> >::iterator itDaughterMap;
[4006893]465 int i, j, code, counter;
466
467 counter = 0;
468 for(i = 0; i < fVertices.size(); ++i)
469 {
470 position = fVertices[i].first;
471 array = fVertices[i].second;
472
473 for(j = 0; j < array->GetEntriesFast(); ++j)
474 {
475 candidate = static_cast<Candidate *>(array->At(j));
476
477 candidate->Position = *position;
478 if(fPositionCoefficient != 1.0)
479 {
480 candidate->Position *= fPositionCoefficient;
481 }
482
[a26a130]483 if(fMomentumCoefficient != 1.0)
484 {
485 candidate->Momentum *= fMomentumCoefficient;
486 }
487
[4006893]488 candidate->M1 = i;
489
490 itDaughterMap = fDaughterMap.find(i);
491 if(itDaughterMap == fDaughterMap.end())
492 {
493 fDaughterMap[i] = make_pair(counter, counter);
494 }
495 else
496 {
497 itDaughterMap->second.second = counter;
498 }
499
500 code = candidate->D1;
501
502 itVertexMap = fInVertexMap.find(code);
503 if(itVertexMap == fInVertexMap.end())
504 {
505 candidate->D1 = -1;
506 }
507 else
508 {
509 code = itVertexMap->second;
510
511 candidate->D1 = code;
512
513 itMotherMap = fMotherMap.find(code);
514 if(itMotherMap == fMotherMap.end())
515 {
516 fMotherMap[code] = make_pair(counter, -1);
517 }
518 else
519 {
520 itMotherMap->second.second = counter;
521 }
522 }
523
524 allParticleOutputArray->Add(candidate);
525
526 ++counter;
[a26a130]527
528 pdgParticle = fPDG->GetParticle(candidate->PID);
529
530 candidate->Charge = pdgParticle ? int(pdgParticle->Charge() / 3.0) : -999;
531
532 if(!pdgParticle) continue;
533
534 pdgCode = TMath::Abs(candidate->PID);
535
536 if(candidate->Status == 1)
537 {
538 stableParticleOutputArray->Add(candidate);
539 }
540 else if(pdgCode <= 5 || pdgCode == 21 || pdgCode == 15)
541 {
542 partonOutputArray->Add(candidate);
543 }
[4006893]544 }
545 }
[d7d2da3]546
547 for(i = 0; i < allParticleOutputArray->GetEntriesFast(); ++i)
548 {
549 candidate = static_cast<Candidate *>(allParticleOutputArray->At(i));
550
[4006893]551 itMotherMap = fMotherMap.find(candidate->M1);
552 if(itMotherMap == fMotherMap.end())
553 {
554 candidate->M1 = -1;
555 candidate->M2 = -1;
556 }
557 else
[d7d2da3]558 {
[4006893]559 candidate->M1 = itMotherMap->second.first;
560 candidate->M2 = itMotherMap->second.second;
[d7d2da3]561 }
[b4786d3]562
[4006893]563 if(candidate->D1 < 0)
[d7d2da3]564 {
565 candidate->D1 = -1;
566 candidate->D2 = -1;
567 }
568 else
569 {
[4006893]570 itDaughterMap = fDaughterMap.find(candidate->D1);
571 if(itDaughterMap == fDaughterMap.end())
572 {
573 candidate->D1 = -1;
574 candidate->D2 = -1;
575 }
576 else
577 {
578 candidate->D1 = itDaughterMap->second.first;
579 candidate->D2 = itDaughterMap->second.second;
580 }
[d7d2da3]581 }
582 }
583}
584
585//---------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.