]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TEvtGen/THepMCParser.cxx
Added HepMC parser. Build it to to a library independent of HepMC. Updated HepMC...
[u/mrichter/AliRoot.git] / TEvtGen / THepMCParser.cxx
1 // module identifier line...
2 // Author: Brian Thorsbro, 24/6-2014
3
4 #include <iostream>
5 #include <sstream>
6 #include <string>
7 #include <list>
8 #include <set>
9 #include <time.h>
10
11 #include "THepMCParser.h"
12 #include "HepMC/IO_GenEvent.h"
13 #include "TObject.h"
14 #include "TTree.h"
15 #include "TClonesArray.h"
16 #include "TFile.h"
17 #include "TParticle.h"
18 #include "TDatabasePDG.h"
19
20 using namespace std;
21
22
23 THepMCParser::THepMCParser(const char * infile) : fTree(0)
24 {
25    HepMC::IO_BaseClass * events = new HepMC::IO_GenEvent(infile, std::ios::in);
26    init(events);
27    delete events;
28    events = 0; // nullptr
29 }
30 THepMCParser::THepMCParser(HepMC::IO_BaseClass * events) : fTree(0)
31 {
32    init(events);
33 }
34 void THepMCParser::init(HepMC::IO_BaseClass * events)
35 {
36    int particlecount = 0;
37    fTree = new TTree("treeEPOS","Tree EPOS");
38    TClonesArray * array = new TClonesArray("TParticle");
39    // array->BypassStreamer();
40    fTree->Branch("Particles",&array); // more flags?
41    THepMCParser::HeavyIonHeader_t heavyIonHeader;
42    fTree->Branch("HeavyIonInfo", &heavyIonHeader, THepMCParser::fgHeavyIonHeaderBranchString);
43    THepMCParser::PdfHeader_t pdfHeader;
44    fTree->Branch("PdfInfo", &pdfHeader, THepMCParser::fgPdfHeaderBranchString);
45    HepMC::GenEvent* evt =  0; // nullptr
46    while ((evt = events->read_next_event())) {
47       string errMsg1 = ParseGenEvent2TCloneArray(evt,array);
48       string errMsg2 = ParseGenEvent2HeaderStructs(evt,heavyIonHeader,pdfHeader);
49       if (errMsg1.length() == 0 && errMsg2.length() == 0) {
50          fTree->Fill();
51       } else {
52          if (errMsg1.length() != 0) cerr << errMsg1 << endl;
53          if (errMsg2.length() != 0) cerr << errMsg2 << endl;
54       }
55       particlecount += array->Capacity();
56    }
57 //   array->Clear();
58 //   delete array;
59 //   array = 0; // nullptr
60    cout << " parsed " << particlecount << " particles" << endl;
61 }
62
63
64 TTree * THepMCParser::GetTTree()
65 {
66    return fTree;
67 }
68 void THepMCParser::WriteTTreeToFile(const char *outfile)
69 {
70    TFile * f = new TFile(outfile, "recreate");
71    fTree->Write();
72    delete f;
73    f = 0; // nullptr
74 }
75
76
77
78 // Based on a validator written by Peter Hristov, CERN
79 bool THepMCParser::IsValidMotherDaughtersConsitency(bool useStdErr, bool requireSecondMotherBeforeDaughters)
80 {
81    bool valid = true;
82    TClonesArray * array = new TClonesArray("TParticle");
83    TBranch* branch = fTree->GetBranch("Particles");
84    branch->SetAddress(&array);
85    Int_t count = branch->GetEntries();
86    for (Int_t idx=0; idx<count; ++idx) {
87       array->Clear();
88       branch->GetEntry(idx); // "fill" the array
89       Int_t nkeep = array->GetEntriesFast();
90       for (Int_t i=0; i<nkeep; i++) {
91          TParticle * part = (TParticle*)array->AddrAt(i);
92          Int_t mum1 = part->GetFirstMother();
93          Int_t mum2 = part->GetSecondMother();
94          Int_t fd = part->GetFirstDaughter();
95          Int_t ld = part->GetLastDaughter();
96          if (mum1>-1 && i<mum1) {
97             valid = false;
98             if (useStdErr) cerr << "Problem: first_mother(" << mum1 << ") after daughter(" << i << ")" << endl;
99          }
100          if (mum2>-1 && i<mum2 && requireSecondMotherBeforeDaughters) {
101             valid = false;
102             if (useStdErr) cerr << "Problem: second_mother(" << mum2 << ") after daughter(" << i << ")" << endl;
103          }
104          if (fd > ld ) {
105             valid = false;
106             if (useStdErr) cerr << "Problem: first_daughter(" << fd << ") > last_daughter(" << ld << ")" << endl;
107          }
108          for (Int_t id=TMath::Max(fd,0); id<=ld; id++) {
109             TParticle * daughter = (TParticle*)array->AddrAt(id);
110             if (daughter->GetFirstMother() != i && daughter->GetSecondMother() != i) {
111                valid = false;
112                if (useStdErr) cerr << "Problem: mother("<< i << ") not registered as mother for either first_daughter("
113                      << daughter->GetFirstMother() << ") or second_daughter("
114                      << daughter->GetSecondMother() << ")" << endl;
115             }
116          }
117       }
118    }
119    delete array;
120    array = 0;
121    return valid;
122 }
123
124 bool THepMCParser::IsValidParticleInvariantMass(bool useStdErr, bool includeStatusCode2Particles)
125 {
126    bool valid = true;
127    TClonesArray *array = new TClonesArray("TParticle");
128    TBranch* branch = fTree->GetBranch("Particles");
129    branch->SetAddress(&array);
130    Int_t count = branch->GetEntries();
131    for (Int_t idx=0; idx<count; ++idx) {
132       array->Clear();
133       branch->GetEntry(idx);
134       Int_t nkeep = array->GetEntries();
135       for (Int_t i=0; i<nkeep; i++) {
136          TParticle * parton = (TParticle*)array->AddrAt(i);
137          if (parton->GetStatusCode()==2 && !includeStatusCode2Particles) {
138             continue;
139          }
140          TLorentzVector v;
141          parton->Momentum(v);
142          Double_t m1 = v.M(); // invariant mass from the particle in the HepMC event
143          TParticlePDG *dbParton = parton->GetPDG();
144          if (!dbParton) {
145             if (useStdErr) cerr << "Warning: could not look up particle with PDG Code: " << parton->GetPdgCode() << endl << endl;
146             continue;
147          }
148          Double_t m2 = dbParton->Mass();
149          bool checkok;
150          if (m2 == 0) {
151             checkok = abs(m1) < 0.0001; // no such thing as negative mass...
152          } else {
153             checkok = abs(1 - m1/m2) < 0.01;
154          }
155          if (!checkok && useStdErr) {
156             cerr << "Problem: " << GetParticleName(parton) << " HepMC:" << m1 << endl;
157             cerr << ListReactionChain(array,i);
158             cerr << endl;
159          }
160          if (!checkok)
161             valid = false;
162       }
163    }
164    delete array;
165    array = 0;
166    return valid;
167 }
168
169 bool THepMCParser::IsValidVertexInvariantMass(bool useStdErr, bool includeStatusCode2Particles)
170 {
171    bool valid = true;
172    TClonesArray * array = new TClonesArray("TParticle");
173    TBranch* branch = fTree->GetBranch("Particles");
174    branch->SetAddress(&array);
175    Int_t count = branch->GetEntries();
176    for (Int_t idx=0; idx<count; ++idx) {
177       array->Clear();
178       branch->GetEntry(idx); // "fill" the array
179       TLorentzVector v_st1;
180       TLorentzVector v_st4;
181       Int_t nkeep = array->GetEntriesFast();
182       for (Int_t i=0; i<nkeep; i++) {
183          TParticle * parton = (TParticle*)array->AddrAt(i);
184          TLorentzVector v_in;
185          parton->Momentum(v_in);
186          if (parton->GetStatusCode()==4) {
187             v_st4 += v_in;
188          } else if (parton->GetStatusCode()==1) {
189             v_st1 += v_in;
190          }
191          if (!includeStatusCode2Particles) { // only check beam particles vs final particles
192             continue;
193          }
194          Int_t fd = parton->GetFirstDaughter();
195          Int_t ld = parton->GetLastDaughter();
196          if (fd == -1) continue; // no daughters, continue loop
197          Int_t mother2 = -1;
198          TLorentzVector v_out;
199          bool oneok = false; bool allok = true; bool agreemother2 = true;
200          ostringstream daughterpdg;
201          ostringstream motherpdg;
202          for (Int_t id=TMath::Max(fd,0); id<=ld; id++) {
203             TParticle * daughter = (TParticle*)array->AddrAt(id);
204             if (fd==id) {
205                daughter->Momentum(v_out);
206                mother2 = daughter->GetSecondMother();
207             } else {
208                TLorentzVector d;
209                daughter->Momentum(d);
210                v_out += d;
211                if (daughter->GetSecondMother() != mother2) agreemother2 = false;
212             }
213             if (daughter->GetFirstMother() == i) {
214                oneok = true;
215             } else {
216                allok = false;
217             }
218             daughterpdg << " " << daughter->GetPdgCode();
219          }
220          motherpdg << " " << parton->GetPdgCode();
221          if (mother2 > -1 && agreemother2) {
222             TParticle * m2 = (TParticle*)array->AddrAt(mother2);
223             TLorentzVector m2v;
224             m2->Momentum(m2v);
225             v_in += m2v;
226             motherpdg << " " << m2->GetPdgCode();
227          }
228          if (oneok && allok && agreemother2) {
229             bool checkok = abs(1 - v_in.M()/v_out.M()) < 0.1;
230             if (!checkok) valid=false;
231             if (!checkok && useStdErr) {
232                //             cerr << "Problem: " << i << "[" << fd << "," << ld << "] PDG:" << motherpdg.str() << " ->" << daughterpdg.str() << "  Inv.mass: " << v_in.M() << " -> " << v_out.M() << endl;
233                cerr << ListReactionChain(array,i);
234                cerr << endl;
235                //             cerr << v_in.Px() << " " << v_in.Py() << " " << v_in.Pz() << " " << v_in.E() << endl;
236             } else if (useStdErr) {
237                //             cerr << "OK: " << i << "[" << fd << "," << ld << "] PDG:" << motherpdg.str() << " ->" << daughterpdg.str() << "  Inv.mass: " << v_in.M() << " -> " << v_out.M() << endl;
238             }
239          }
240       }
241       bool checkok = abs(1 - v_st4.M()/v_st1.M()) < 0.001;
242       if (!checkok) valid=false;
243       if (!checkok && useStdErr) {
244          cerr << " BEAM PARTICLES -> FINAL PARTICLES " << endl;
245          cerr << " status 4 (" << v_st4.M() << ") -> status 1 (" << v_st1.M() << ")" << endl << endl;
246       }
247    }
248    delete array;
249    array = 0;
250    return valid;
251 }
252
253 string THepMCParser::GetParticleName(TParticle * thisPart)
254 {
255    TParticlePDG *dbPart = thisPart->GetPDG();
256    ostringstream name;
257    if (dbPart) {
258       name << dbPart->GetName() << "(" << dbPart->PdgCode() << "," << dbPart->Mass() << ")";
259    } else {
260       name << thisPart->GetPdgCode() << "(NoDBinfo)";
261    }
262    return name.str();
263 }
264
265 string THepMCParser::ListReactionChain(TClonesArray * particles, Int_t particleId)
266 {
267    ostringstream output;
268
269    TParticle * part = (TParticle*)particles->AddrAt(particleId);
270    Int_t m1id = part->GetFirstMother();
271    Int_t m2id = part->GetSecondMother();
272    if (m1id > 1) { // ignore the initial collision with beam particles
273       ostringstream inStr;
274       ostringstream outStr;
275       TParticle * m1 = (TParticle*)particles->AddrAt(m1id);
276       TLorentzVector v_in;
277       m1->Momentum(v_in);
278       inStr << GetParticleName(m1) << "[" << v_in.M() << "]";
279       if (m2id > 1) {
280          TParticle * m2 = (TParticle*)particles->AddrAt(m2id);
281          TLorentzVector v_m2;
282          m2->Momentum(v_m2);
283          v_in += v_m2;
284          inStr << " " << GetParticleName(m2) << "[" << v_m2.M() << "]";
285       }
286       Int_t fd = m1->GetFirstDaughter();
287       Int_t ld = m1->GetLastDaughter();
288       TLorentzVector v_out;
289       part->Momentum(v_out);
290       outStr << GetParticleName(part) << "[" << v_out.M() << "]";
291       for (Int_t i=fd; i<=ld; ++i) {
292          if (i!=particleId) {
293             TParticle * d = (TParticle*)particles->AddrAt(i);
294             TLorentzVector v_d;
295             d->Momentum(v_d);
296             v_out += v_d;
297             outStr << " " << GetParticleName(d) << "[" << v_d.M() << "]";
298          }
299       }
300       output << "Parent reaction, inv mass: " << v_in.M() << " -> " << v_out.M() << endl;
301       output << " - partons: " << inStr.str() << " -> " << outStr.str() << endl;
302    }
303    Int_t fd = part->GetFirstDaughter();
304    Int_t ld = part->GetLastDaughter();
305    if (fd > -1) {
306       ostringstream inStr;
307       ostringstream outStr;
308       TLorentzVector v_in;
309       part->Momentum(v_in);
310       inStr << GetParticleName(part) << "[" << v_in.M() << "]";
311
312       TParticle * f = (TParticle*)particles->AddrAt(fd);
313       m2id = f->GetSecondMother();
314       if (m2id == particleId) {
315          m2id = f->GetFirstMother();
316       }
317       if (m2id > -1) {
318          TParticle * m2 = (TParticle*)particles->AddrAt(m2id);
319          TLorentzVector v_m2;
320          m2->Momentum(v_m2);
321          v_in += v_m2;
322          inStr << " " << GetParticleName(m2) << "[" << v_m2.M() << "]";
323       }
324       TLorentzVector v_out;
325       f->Momentum(v_out);
326       outStr << GetParticleName(f) << "[" << v_out.M() << "]";
327       for (Int_t i=fd+1; i<=ld; ++i) {
328          TParticle * d = (TParticle*)particles->AddrAt(i);
329          TLorentzVector v_d;
330          d->Momentum(v_d);
331          v_out += v_d;
332          outStr << " " << GetParticleName(d) << "[" << v_d.M() << "]";
333       }
334       output << "Child reaction, inv mass: " << v_in.M() << " -> " << v_out.M() << endl;
335       output << " - partons: " << inStr.str() << " -> " << outStr.str() << endl;
336    } else {
337       output << "Child reaction" << endl << " - none" << endl;
338    }
339
340    return output.str();
341 }
342
343
344 string THepMCParser::ParseGenEvent2TCloneArray(HepMC::GenEvent * genEvent, TClonesArray * array, bool requireSecondMotherBeforeDaughters)
345 {
346    if (requireSecondMotherBeforeDaughters) {
347       return "requireSecondMotherBeforeDaughters not implemented yet!";
348    }
349    array->Clear();
350    ostringstream errMsgStream;
351    map<int,Int_t> partonMemory; // unordered_map in c++11 - but probably not much performance gain from that: log(n) vs log(1) where constant can be high
352
353    // Check event with HepMC's internal validation algorithm
354    if (!genEvent->is_valid()) {
355       errMsgStream << "Error with event id: " << genEvent->event_number() << ", event is not valid!";
356       return errMsgStream.str();
357    }
358
359    // Pull out the beam particles from the event
360    const pair<HepMC::GenParticle *,HepMC::GenParticle *> beamparts = genEvent->beam_particles();
361
362    // Four sanity checks:
363    // - Beam particles exists and are not the same
364    // - Both beam particles should have no production vertices, they come from the beams
365    // - Both beam particles should have defined end vertices, as they both should contribute
366    // - Both beam particles should have the exact same end vertex
367    if (!beamparts.first || !beamparts.second || beamparts.first->barcode()==beamparts.second->barcode()) {
368       errMsgStream << "Error with event id: " << genEvent->event_number() << ", beam particles doesn't exists or are the same";
369       return errMsgStream.str();
370    }
371    if (beamparts.first->production_vertex() || beamparts.second->production_vertex()) {
372       errMsgStream << "Error with event id: " << genEvent->event_number() << ", beam particles have production vertex/vertices...";
373       return errMsgStream.str();
374    }
375    if (!beamparts.first->end_vertex() || !beamparts.second->end_vertex()) {
376       errMsgStream << "Error with event id: " << genEvent->event_number() << ", beam particles have undefined end vertex/vertices...";
377       return errMsgStream.str();
378    }
379    if (beamparts.first->end_vertex() != beamparts.second->end_vertex()) {
380       errMsgStream << "Error with event id: " << genEvent->event_number() << ", beam particles do not collide in the same end vertex.";
381       return errMsgStream.str();
382    }
383
384    // Set the array to hold the number of particles in the event
385    array->Expand(genEvent->particles_size());
386
387    // Create a TParticle for each beam particle
388    new((*array)[0]) TParticle(
389          beamparts.first->pdg_id(),
390          beamparts.first->status(), // check if status has the same meaning
391          -1, // no mother1
392          -1, // no mother2
393          -1, // first daughter not known yet
394          -1, // last daughter not known yet
395          beamparts.first->momentum().px(),
396          beamparts.first->momentum().py(),
397          beamparts.first->momentum().pz(),
398          beamparts.first->momentum().e(),
399          0, // no production vertex, so zero?
400          0,
401          0,
402          0
403    );
404    partonMemory[beamparts.first->barcode()] = 0;
405    new((*array)[1]) TParticle(
406          beamparts.second->pdg_id(),
407          beamparts.second->status(),
408          -1, // no mother1
409          -1, // no mother2
410          -1, // first daughter not known yet
411          -1, // last daughter not known yet
412          beamparts.second->momentum().px(),
413          beamparts.second->momentum().py(),
414          beamparts.second->momentum().pz(),
415          beamparts.second->momentum().e(),
416          0, // no production vertex, so zero?
417          0,
418          0,
419          0
420    );
421    partonMemory[beamparts.second->barcode()] = 1;
422    Int_t arrayID = 2; // start counting IDs after the beam particles
423
424    // Do first vertex as a special case for performance, since its often has the most daughters and both mothers are known
425    Int_t firstDaughter = arrayID;
426    for (HepMC::GenVertex::particles_out_const_iterator iter = beamparts.first->end_vertex()->particles_out_const_begin();
427          iter != beamparts.first->end_vertex()->particles_out_const_end();
428          ++iter) {
429       new((*array)[arrayID]) TParticle(
430             (*iter)->pdg_id(),
431             (*iter)->status(),
432             0, // beam particle 1
433             1, // beam particle 2
434             -1, // first daughter not known yet
435             -1, // last daughter not known yet
436             (*iter)->momentum().px(),
437             (*iter)->momentum().py(),
438             (*iter)->momentum().pz(),
439             (*iter)->momentum().e(),
440             beamparts.first->end_vertex()->position().x(),
441             beamparts.first->end_vertex()->position().y(),
442             beamparts.first->end_vertex()->position().z(),
443             beamparts.first->end_vertex()->position().t()
444       );
445       partonMemory[(*iter)->barcode()] = arrayID;
446       ++arrayID;
447    }
448    Int_t lastDaughter = arrayID-1;
449    ((TParticle*)array->AddrAt(0))->SetFirstDaughter(firstDaughter); // beam particle 1
450    ((TParticle*)array->AddrAt(0))->SetLastDaughter(lastDaughter);
451    ((TParticle*)array->AddrAt(1))->SetFirstDaughter(firstDaughter); // beam particle 2
452    ((TParticle*)array->AddrAt(1))->SetLastDaughter(lastDaughter);
453
454    // Build vertex list by exploring tree and sorting such that daughters comes after mothers
455    list<HepMC::GenVertex*> vertexList;
456    set<int> vertexSearchSet;
457    ExploreVertex(beamparts.first->end_vertex(),vertexList,vertexSearchSet,requireSecondMotherBeforeDaughters);
458
459    // Analyze each vertex
460    for (list<HepMC::GenVertex*>::iterator i = vertexList.begin(); i != vertexList.end(); ++i) {
461       HepMC::GenVertex * vertex = (*i);
462
463       // first establish mother relations
464       HepMC::GenVertex::particles_in_const_iterator iter = vertex->particles_in_const_begin();
465       if (iter == vertex->particles_in_const_end()) {
466          return "Particle without a mother, and its not a beam particle!";
467       }
468       int motherA = partonMemory[(*iter)->barcode()];
469       if (((TParticle*)array->AddrAt(motherA))->GetFirstDaughter() > -1) {
470          return "Trying to assign new daughters to a particle that already has daughters defined!";
471       }
472       ++iter;
473       int motherB = -1;
474       if (iter != vertex->particles_in_const_end()) {
475          motherB = partonMemory[(*iter)->barcode()];
476          if (((TParticle*)array->AddrAt(motherB))->GetFirstDaughter() > -1) {
477             return "Trying to assign new daughters to a particle that already has daughters defined!";
478          }
479          ++iter;
480          if (iter != vertex->particles_in_const_end()) {
481             return "Particle with more than two mothers!";
482          }
483       }
484       if (motherB > -1 && motherB < motherA) {
485          int swap = motherA; motherA = motherB; motherB = swap;
486       }
487
488       // add the particles to the array, important that they are add in succession with respect to arrayID
489       firstDaughter = arrayID;
490       for (iter = vertex->particles_out_const_begin();
491             iter != vertex->particles_out_const_end();
492             ++iter) {
493          new((*array)[arrayID]) TParticle(
494                (*iter)->pdg_id(),
495                (*iter)->status(),
496                motherA, // mother 1
497                motherB, // mother 2
498                -1, // first daughter, if applicable, not known yet
499                -1, // last daughter, if applicable, not known yet
500                (*iter)->momentum().px(),
501                (*iter)->momentum().py(),
502                (*iter)->momentum().pz(),
503                (*iter)->momentum().e(),
504                vertex->position().x(),
505                vertex->position().y(),
506                vertex->position().z(),
507                vertex->position().t()
508          );
509          partonMemory[(*iter)->barcode()] = arrayID;
510          ++arrayID;
511       }
512       lastDaughter = arrayID-1;
513       if (lastDaughter < firstDaughter) {
514          return "Vertex with no out particles, should not be possible!";
515       }
516       // update mother with daughter interval
517       ((TParticle*)array->AddrAt(motherA))->SetFirstDaughter(firstDaughter);
518       ((TParticle*)array->AddrAt(motherA))->SetLastDaughter(lastDaughter);
519       if (motherB > -1) {
520          ((TParticle*)array->AddrAt(motherB))->SetFirstDaughter(firstDaughter);
521          ((TParticle*)array->AddrAt(motherB))->SetLastDaughter(lastDaughter);
522       }
523    }
524
525    return "";
526 }
527
528
529 void THepMCParser::ExploreVertex(HepMC::GenVertex * vertex, list<HepMC::GenVertex*> & vertexList, set<int> & vertexSearchSet, bool requireSecondMotherBeforeDaughters)
530 {
531    for (HepMC::GenVertex::particles_out_const_iterator iter = vertex->particles_out_const_begin();
532          iter != vertex->particles_out_const_end();
533          ++iter) {
534       HepMC::GenVertex * testVertex = (*iter)->end_vertex();
535       if (testVertex) {
536          bool foundVertex = vertexSearchSet.find((*iter)->end_vertex()->barcode()) != vertexSearchSet.end();
537          if (requireSecondMotherBeforeDaughters) {
538             // redo this algorithem to move subtree instead of node....
539             // its not completely error proof in its current implementation even though the error is extremely rare
540
541             if (foundVertex) for (list<HepMC::GenVertex*>::iterator i = vertexList.begin(); i != vertexList.end(); ++i) {
542                if ((*i)->barcode() == testVertex->barcode()) {
543                   vertexList.erase(i);
544                   cout << " it happened, the vertex parsing order had to be changed " << endl;
545                   break;
546                }
547             } else {
548                vertexSearchSet.insert((*iter)->end_vertex()->barcode());
549             }
550             vertexList.push_back(testVertex);
551             if (!foundVertex) ExploreVertex(testVertex,vertexList,vertexSearchSet,requireSecondMotherBeforeDaughters);
552
553          } else {
554             if (!foundVertex) {
555                vertexSearchSet.insert((*iter)->end_vertex()->barcode());
556                vertexList.push_back(testVertex);
557                ExploreVertex(testVertex,vertexList,vertexSearchSet,requireSecondMotherBeforeDaughters);
558             }
559          }
560       }
561    }
562 }
563
564
565
566 const char * THepMCParser::fgHeavyIonHeaderBranchString = "Ncoll_hard/I:Npart_proj:Npart_targ:Ncoll:spectator_neutrons:spectator_protons:N_Nwounded_collisions:Nwounded_N_collisions:Nwounded_Nwounded_collisions:impact_parameter/F:event_plane_angle:eccentricity:sigma_inel_NN";
567 const char * THepMCParser::fgPdfHeaderBranchString = "id1/I:id2:pdf_id1:pdf_id2:x1/D:x2:scalePDF:pdf1:pdf2";
568
569 string THepMCParser::ParseGenEvent2HeaderStructs(HepMC::GenEvent * genEvent, HeavyIonHeader_t & heavyIonHeader, PdfHeader_t & pdfHeader, bool fillZeroOnMissingHeavyIon, bool fillZeroOnMissingPdf)
570 {
571    HepMC::HeavyIon * heavyIon = genEvent->heavy_ion();
572    HepMC::PdfInfo * pdfInfo = genEvent->pdf_info();
573    if ((!heavyIon && !fillZeroOnMissingHeavyIon) || (!pdfInfo && !fillZeroOnMissingPdf)) {
574       return "HeavyIonInfo and/or PdfInfo not defined for this event, did you read it with IO_GenEvent?";
575    }
576    if (heavyIon) {
577       heavyIonHeader.Ncoll_hard = heavyIon->Ncoll_hard();
578       heavyIonHeader.Npart_proj = heavyIon->Npart_proj();
579       heavyIonHeader.Npart_targ = heavyIon->Npart_targ();
580       heavyIonHeader.Ncoll = heavyIon->Ncoll();
581       heavyIonHeader.spectator_neutrons = heavyIon->spectator_neutrons();
582       heavyIonHeader.spectator_protons = heavyIon->spectator_protons();
583       heavyIonHeader.N_Nwounded_collisions = heavyIon->N_Nwounded_collisions();
584       heavyIonHeader.Nwounded_N_collisions = heavyIon->Nwounded_N_collisions();
585       heavyIonHeader.Nwounded_Nwounded_collisions = heavyIon->Nwounded_Nwounded_collisions();
586       heavyIonHeader.impact_parameter = heavyIon->impact_parameter();
587       heavyIonHeader.event_plane_angle = heavyIon->event_plane_angle();
588       heavyIonHeader.eccentricity = heavyIon->eccentricity();
589       heavyIonHeader.sigma_inel_NN = heavyIon->sigma_inel_NN();
590    } else {
591       heavyIonHeader.Ncoll_hard = 0;
592       heavyIonHeader.Npart_proj = 0;
593       heavyIonHeader.Npart_targ = 0;
594       heavyIonHeader.Ncoll = 0;
595       heavyIonHeader.spectator_neutrons = 0;
596       heavyIonHeader.spectator_protons = 0;
597       heavyIonHeader.N_Nwounded_collisions = 0;
598       heavyIonHeader.Nwounded_N_collisions = 0;
599       heavyIonHeader.Nwounded_Nwounded_collisions = 0;
600       heavyIonHeader.impact_parameter = 0.0;
601       heavyIonHeader.event_plane_angle = 0.0;
602       heavyIonHeader.eccentricity = 0.0;
603       heavyIonHeader.sigma_inel_NN = 0.0;
604    }
605    if (pdfInfo) {
606       pdfHeader.id1 = pdfInfo->id1();
607       pdfHeader.id2 = pdfInfo->id2();
608       pdfHeader.pdf_id1 = pdfInfo->pdf_id1();
609       pdfHeader.pdf_id2 = pdfInfo->pdf_id2();
610       pdfHeader.x1 = pdfInfo->x1();
611       pdfHeader.x2 = pdfInfo->x2();
612       pdfHeader.scalePDF = pdfInfo->scalePDF();
613       pdfHeader.pdf1 = pdfInfo->pdf1();
614       pdfHeader.pdf2 = pdfInfo->pdf2();
615    } else {
616       pdfHeader.id1 = 0;
617       pdfHeader.id2 = 0;
618       pdfHeader.pdf_id1 = 0;
619       pdfHeader.pdf_id2 = 0;
620       pdfHeader.x1 = 0.0;
621       pdfHeader.x2 = 0.0;
622       pdfHeader.scalePDF = 0.0;
623       pdfHeader.pdf1 = 0.0;
624       pdfHeader.pdf2 = 0.0;
625    }
626    return "";
627 }
628
629
630
631
632
633
634