15 #include <TLorentzVector.h>
18 #include <TRootIOCtor.h>
20 #include "Framework/Conventions/GBuild.h"
34 using std::setprecision;
38 using namespace genie;
54 TClonesArray(
"genie::GHepParticle")
60 TClonesArray(
"genie::GHepParticle", size)
66 TClonesArray(
"genie::GHepParticle", record.GetEntries())
73 TClonesArray(
"genie::GHepParticle"),
94 LOG(
"GHEP",
pWARN) <<
"Returning NULL interaction";
108 if( position >=0 && position < this->GetEntries() ) {
110 if(particle)
return particle;
113 <<
"No particle found with: (pos = " << position <<
")";
124 int nentries = this->GetEntries();
125 for(
int i = start; i < nentries; i++) {
127 if(p->
Status() == status && p->
Pdg() == pdg)
return p;
131 <<
"No particle found with: (pos >= " << start
132 <<
", pdg = " << pdg <<
", ist = " << status <<
")";
143 int nentries = this->GetEntries();
144 for(
int i = start; i < nentries; i++) {
146 if(p->
Status() == status && p->
Pdg() == pdg)
return i;
150 <<
"No particle found with: (pos >= " << start
151 <<
", pdg = " << pdg <<
", ist = " << status <<
")";
161 int nentries = this->GetEntries();
162 for(
int i = start; i < nentries; i++) {
164 if( p->
Compare(particle) )
return i;
168 <<
"No particle found with pos >= " << start
169 <<
" matching particle: " << *particle;
180 int nentries = this->GetEntries();
181 if(position<0 || position>=nentries)
return 0;
183 vector<int> * descendants =
new vector<int>;
187 descendants->push_back(position);
191 for(
int i = 0; i < nentries; i++) {
192 if(i==position)
continue;
195 bool is_descendant=
false;
198 if(mom==position) is_descendant=
true;
200 descendants->push_back(i);
216 int p0pdg = p0->
Pdg();
218 int p1pdg = p1->
Pdg();
289 if(ipos>-1)
return this->
Particle(ipos);
299 if(ipos>-1)
return this->
Particle(ipos);
309 if(ipos>-1)
return this->
Particle(ipos);
319 if(ipos>-1)
return this->
Particle(ipos);
329 if(ipos>-1)
return this->
Particle(ipos);
338 if(ipos>-1)
return this->
Particle(ipos);
348 if(ipos>-1)
return this->
Particle(ipos);
408 if(dau1==-1 && dau2==-1)
return -1;
410 for(
int i=dau1; i<=dau2; i++) {
412 int dpdgc = dp->
Pdg();
428 int ipos = (nucleus) ? 2 : 1;
436 if(isN && p->
Status()==ist)
return ipos;
448 int ipos = (nucleus) ? 2 : 1;
465 if(!probe)
return -1;
479 unsigned int nentries = 0;
481 for(
int i = start; i < this->GetEntries(); i++) {
483 if(p->
Pdg()==pdg && p->
Status()==ist) nentries++;
490 unsigned int nentries = 0;
492 for(
int i = start; i < this->GetEntries(); i++) {
494 if(p->
Pdg()==pdg) nentries++;
503 unsigned int pos = this->GetEntries();
505 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
507 <<
"Adding particle with pdgc = " << p.
Pdg() <<
" at slot = " << pos;
517 int pdg,
GHepStatus_t status,
int mom1,
int mom2,
int dau1,
int dau2,
518 const TLorentzVector & p,
const TLorentzVector & v)
522 unsigned int pos = this->GetEntries();
524 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
526 <<
"Adding particle with pdgc = " << pdg <<
" at slot = " << pos;
528 new ((*this)[pos])
GHepParticle(pdg,status, mom1,mom2,dau1,dau2, p, v);
536 int pdg,
GHepStatus_t status,
int mom1,
int mom2,
int dau1,
int dau2,
537 double px,
double py,
double pz,
double E,
538 double x,
double y,
double z,
double t)
542 unsigned int pos = this->GetEntries();
544 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
546 <<
"Adding particle with pdgc = " << pdg <<
" at slot = " << pos;
549 pdg, status, mom1, mom2, dau1, dau2, px, py, pz, E, x, y, z, t);
558 int pos = this->GetEntries() - 1;
560 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
562 <<
"Updating the daughter-list for the mother of particle at: " << pos;
569 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
570 LOG(
"GHEP",
pINFO) <<
"Mother particle is at slot: " << mom_pos;
572 if(mom_pos==-1)
return;
583 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
585 <<
"Done! Daughter-list is compact: [" << pos <<
", " << pos <<
"]";
593 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
595 <<
"Done! Daughter-list is compact: [" << pos <<
", " << dau2 <<
"]";
603 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
605 <<
"Done! Daughter-list is compact: [" << dau1 <<
", " << pos <<
"]";
613 <<
"Daughter-list is not compact - Running compactifier";
619 LOG(
"GHEP",
pNOTICE) <<
"Removing all intermediate particles from GHEP";
640 <<
"Removing: " << p->
Name() <<
" from slot: " << i;
645 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
646 LOG(
"GHEP",
pDEBUG) <<
"Compressing GHEP record to remove empty slots";
653 int n = this->GetEntries();
662 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
664 <<
"Particle's " << i <<
" daughter list is "
665 << ((compact) ?
"compact" :
"__not__ compact");
673 int ndau = dau2-dau1+1;
675 if(dau1==-1) {ndau=0;}
678 while (curr_pos > ndp) {
693 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
695 <<
"Done ompactifying daughter-list for particle at slot: " << i;
702 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
703 LOG(
"GHEP",
pDEBUG) <<
"Examining daughter-list of particle at: " << pos;
705 vector<int> daughters;
712 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
713 LOG(
"GHEP",
pDEBUG) <<
"Particle at: " << i <<
" is a daughter";
715 daughters.push_back(i);
720 bool is_compact =
true;
721 if(daughters.size()>1) {
722 sort(daughters.begin(), daughters.end());
723 vector<int>::iterator diter = daughters.begin();
725 for(; diter != daughters.end(); ++diter) {
727 is_compact = is_compact && (TMath::Abs(prev-curr)<=1);
731 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
733 <<
"Daughter-list of particle at: " << pos <<
" is "
734 << (is_compact ?
"" :
"not") <<
" compact";
754 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
755 LOG(
"GHEP",
pINFO) <<
"Swapping GHepParticles : " << i <<
" <--> " << j;
757 int n = this->GetEntries();
758 assert(i>=0 && j>=0 && i<n && j<n);
773 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
775 << pi->
Name() <<
"(previously at pos: " << j
776 <<
") is now at pos: " << i <<
" -> Notify daughters";
778 for(
int k=0; k<n; k++) {
779 if(this->
Particle(k)->FirstMother()==j) {
786 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
788 << pj->
Name() <<
"(previously at pos: " << i
789 <<
") is now at pos: " << j <<
" -> Notify daughters";
791 for(
int k=0; k<n; k++) {
792 if(this->
Particle(k)->FirstMother()==i) {
816 dau1 = (dau1<0) ? i2 : TMath::Min(dau1,i2);
817 dau2 = (dau2<0) ? i2 : TMath::Max(dau2,i2);
822 p1 -> SetFirstDaughter (dau1);
823 p1 -> SetLastDaughter (dau2);
829 fVtx->SetXYZT(x,y,z,t);
834 fVtx->SetXYZT(vtx.X(),vtx.Y(),vtx.Z(),vtx.T());
839 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
840 LOG(
"GHEP",
pDEBUG) <<
"Initializing GHepRecord";
848 fVtx =
new TLorentzVector(0,0,0,0);
863 this->SetOwner(
true);
868 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
869 LOG(
"GHEP",
pDEBUG) <<
"Cleaning up GHepRecord";
876 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
877 LOG(
"GHEP",
pDEBUG) <<
"Reseting GHepRecord";
897 TClonesArray::Clear(opt);
913 unsigned int ientry = 0;
915 TIter ghepiter(&record);
927 TLorentzVector * v = record.
Vertex();
928 fVtx->SetXYZT(v->X(),v->Y(),v->Z(),v->T());
951 TBits bitwiseand = flags & mask;
952 bool accept = (bitwiseand.CountBits() == 0);
977 bool accept_input_print_level =
981 int printlevel = (accept_input_print_level) ?
fPrintLevel : 3;
982 int printlevel_orig = printlevel;
984 bool showpos =
false;
985 if(printlevel >= 10) {
991 stream << setfill(
'-') << setw(115) <<
"|";
993 stream <<
"\n|GENIE GHEP Event Record [print level: "
994 << setfill(
' ') << setw(3) << printlevel_orig <<
"]"
995 << setfill(
' ') << setw(73) <<
"|";
998 stream << setfill(
'-') << setw(115) <<
"|";
1001 stream << setfill(
' ') << setw(6) <<
"Idx | "
1002 << setfill(
' ') << setw(16) <<
"Name | "
1003 << setfill(
' ') << setw(6) <<
"Ist | "
1004 << setfill(
' ') << setw(13) <<
"PDG | "
1005 << setfill(
' ') << setw(12) <<
"Mother | "
1006 << setfill(
' ') << setw(12) <<
"Daughter | "
1007 << setfill(
' ') << setw(10) << ((showpos) ?
"Px(x) |" :
"Px | ")
1008 << setfill(
' ') << setw(10) << ((showpos) ?
"Py(y) |" :
"Py | ")
1009 << setfill(
' ') << setw(10) << ((showpos) ?
"Pz(z) |" :
"Pz | ")
1010 << setfill(
' ') << setw(10) << ((showpos) ?
"E(t) |" :
"E | ")
1011 << setfill(
' ') << setw(10) <<
"m | ";
1014 stream << setfill(
'-') << setw(115) <<
"|";
1017 TObjArrayIter piter(
this);
1018 TVector3 polarization(0,0,0);
1020 unsigned int idx = 0;
1030 stream << setfill(
' ') << setw(3) << idx++ <<
" | ";
1031 stream << setfill(
' ') << setw(13) << p->
Name() <<
" | ";
1032 stream << setfill(
' ') << setw(3) << p->
Status() <<
" | ";
1033 stream << setfill(
' ') << setw(10) << p->
Pdg() <<
" | ";
1034 stream << setfill(
' ') << setw(3) << p->
FirstMother() <<
" | ";
1035 stream << setfill(
' ') << setw(3) << p->
LastMother() <<
" | ";
1036 stream << setfill(
' ') << setw(3) << p->
FirstDaughter() <<
" | ";
1037 stream << setfill(
' ') << setw(3) << p->
LastDaughter() <<
" | ";
1038 stream << std::fixed << setprecision(3);
1039 stream << setfill(
' ') << setw(7) << p->
Px() <<
" | ";
1040 stream << setfill(
' ') << setw(7) << p->
Py() <<
" | ";
1041 stream << setfill(
' ') << setw(7) << p->
Pz() <<
" | ";
1042 stream << setfill(
' ') << setw(7) << p->
E() <<
" | ";
1045 stream << setfill(
' ') << setw(7) << p->
Mass() <<
" | ";
1047 stream << setfill(
'*') << setw(7) << p->
Mass() <<
" | M = "
1048 << p->
P4()->M() <<
" ";
1052 stream <<
"P = (" << polarization.x() <<
"," << polarization.y()
1053 <<
"," << polarization.z() <<
")";
1063 stream << setfill(
' ') << setw(6) <<
" | ";
1064 stream << setfill(
' ') << setw(16) <<
" | ";
1065 stream << setfill(
' ') << setw(6) <<
" | ";
1066 stream << setfill(
' ') << setw(13) <<
" | ";
1067 stream << setfill(
' ') << setw(6) <<
" | ";
1068 stream << setfill(
' ') << setw(6) <<
" | ";
1069 stream << setfill(
' ') << setw(6) <<
" | ";
1070 stream << setfill(
' ') << setw(6) <<
" | ";
1071 stream << std::fixed << setprecision(3);
1072 stream << setfill(
' ') << setw(7) << p->
Vx() <<
" | ";
1073 stream << setfill(
' ') << setw(7) << p->
Vy() <<
" | ";
1074 stream << setfill(
' ') << setw(7) << p->
Vz() <<
" | ";
1075 stream << setfill(
' ') << setw(7) << p->
Vt() <<
" | ";
1076 stream << setfill(
' ') << setw(10) <<
" | ";
1102 stream << setfill(
'-') << setw(115) <<
"|";
1106 stream << setfill(
' ') << setw(17) <<
"Fin-Init: "
1107 << setfill(
' ') << setw(6) <<
" "
1108 << setfill(
' ') << setw(18) <<
" "
1109 << setfill(
' ') << setw(12) <<
" "
1110 << setfill(
' ') << setw(12) <<
" | ";
1111 stream << std::fixed << setprecision(3);
1112 stream << setfill(
' ') << setw(7) << sum_px <<
" | ";
1113 stream << setfill(
' ') << setw(7) << sum_py <<
" | ";
1114 stream << setfill(
' ') << setw(7) << sum_pz <<
" | ";
1115 stream << setfill(
' ') << setw(7) << sum_E <<
" | ";
1116 stream << setfill(
' ') << setw(10) <<
" | ";
1119 stream << setfill(
'-') << setw(115) <<
"|";
1126 stream << setfill(
' ') << setw(17) <<
"Vertex: ";
1127 stream << setfill(
' ') << setw(11)
1128 << ((probe) ? probe->
Name() :
"unknown probe") <<
" @ (";
1130 stream << std::fixed << setprecision(5);
1131 stream <<
"x = " << setfill(
' ') << setw(11) << this->
Vertex()->X() <<
" m, ";
1132 stream <<
"y = " << setfill(
' ') << setw(11) << this->
Vertex()->Y() <<
" m, ";
1133 stream <<
"z = " << setfill(
' ') << setw(11) << this->
Vertex()->Z() <<
" m, ";
1134 stream << std::scientific << setprecision(6);
1135 stream <<
"t = " << setfill(
' ') << setw(15) << this->
Vertex()->T() <<
" s) ";
1136 stream << std::fixed << setprecision(3);
1137 stream << setfill(
' ') << setw(2) <<
"|";
1140 stream << setfill(
'-') << setw(115) <<
"|";
1147 stream <<
"Err flag [bits:" <<
fEventFlags->GetNbits()-1 <<
"->0] : "
1149 <<
"1st set: " << setfill(
' ') << setw(56)
1154 stream <<
"Err mask [bits:" <<
fEventMask->GetNbits()-1 <<
"->0] : "
1156 <<
"Is unphysical: " << setfill(
' ') << setw(5)
1158 <<
"Accepted: " << setfill(
' ') << setw(5)
1162 stream << setfill(
'-') << setw(115) <<
"|";
1167 stream << std::scientific << setprecision(5);
1169 stream <<
"sig(Ev) = "
1175 stream <<
" dsig(y;E)/dy = " << setfill(
' ') << setw(13) <<
fDiffXSec/
units::cm2 <<
" cm^2 |";
1178 stream <<
" d2sig(x,y;E)/dxdy = " << setfill(
' ') << setw(13) <<
fDiffXSec/
units::cm2 <<
" cm^2 |";
1181 stream <<
" d3sig(x,y,t;E)/dxdydt = " << setfill(
' ') << setw(13) <<
fDiffXSec/
units::cm2 <<
" cm^2/GeV^2 |";
1184 stream <<
" dsig(Q2;E)/dQ2 = " << setfill(
' ') << setw(13) <<
fDiffXSec/
units::cm2 <<
" cm^2/GeV^2 |";
1187 stream <<
" dsig(Q2,v,p;E)/dQ2dvdp =" << setfill(
' ') << setw(13) <<
fDiffXSec/
units::cm2 <<
" cm^2/GeV^4 |";
1190 stream <<
" d2sig(W,Q2;E)/dWdQ2 = " << setfill(
' ') << setw(13) <<
fDiffXSec/
units::cm2 <<
" cm^2/GeV^3 |";
1193 stream <<
" d4sig(W,Q2,cos(theta),phi;E)/dWdQ2dOmega = " << setfill(
' ') << setw(13) <<
fDiffXSec/
units::cm2 <<
" cm^2/GeV^3 |";
1196 stream <<
" dsig(Ev;{K_s})/dK = " << setfill(
' ') << setw(13) <<
fDiffXSec/
units::cm2 <<
" cm^2/{K} |";
1198 stream <<
" Weight = "
1199 << setfill(
' ') << setw(16)
1200 << std::fixed << setprecision(5)
1205 stream << setfill(
'-') << setw(115) <<
"|";
1209 stream << setfill(
' ');
1213 else stream <<
"NULL Interaction!" << endl;
static void SetPrintLevel(int print_level)
void SetFirstMother(int m)
int RescatterCode(void) const
virtual GHepParticle * Particle(int position) const
virtual int ParticlePosition(int pdg, GHepStatus_t i, int start=0) const
virtual void Copy(const GHepRecord &record)
virtual void UpdateDaughterLists(void)
double fProb
event probability (for given flux neutrino and density-weighted path-length for target element) ...
void Print(ostream &stream) const
bool IsOnMassShell(void) const
virtual void SwapParticles(int i, int j)
double E(void) const
Get energy.
virtual GHepParticle * HitElectron(void) const
virtual Interaction * Summary(void) const
virtual GHepParticle * FindParticle(int pdg, GHepStatus_t ist, int start) const
const TLorentzVector * P4(void) const
int FirstDaughter(void) const
virtual int RemnantNucleusPosition(void) const
string BoolAsYNString(bool b)
enum genie::EGHepFlag GHepFlag_t
virtual void RemoveIntermediateParticles(void)
virtual void FinalizeDaughterLists(void)
TLorentzVector * fVtx
vertex in the detector coordinate system
virtual int HitNucleonPosition(void) const
bool IsDarkMatter(int pdgc)
TBits * fEventMask
an input bit-field mask allowing one to ignore bits set in fEventFlags
virtual unsigned int NEntries(int pdg, GHepStatus_t ist, int start=0) const
double Mass(void) const
Mass that corresponds to the PDG code.
double Pz(void) const
Get Pz.
static const char * Describe(GHepFlag_t flag)
GHepStatus_t Status(void) const
virtual void AttachSummary(Interaction *interaction)
virtual TLorentzVector * Vertex(void) const
double Px(void) const
Get Px.
bool Compare(const GHepParticle *p) const
virtual GHepParticle * Probe(void) const
int LastMother(void) const
double Vt(void) const
Get production time.
int FirstMother(void) const
string Name(void) const
Name that corresponds to the PDG code.
Enumeration of GENIE event generation modes.
double fXSec
cross section for selected event
Summary information for an interaction.
Interaction * fInteraction
attached summary information
int LastDaughter(void) const
TBits * fEventFlags
event flags indicating various pathologies or an unphysical event
static unsigned int NFlags(void)
virtual bool Accept(void) const
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
void SetLastDaughter(int d)
bool IsAntiDarkMatter(int pdgc)
static constexpr double cm2
bool HasDaughters(void) const
virtual bool HasCompactDaughterList(int pos)
void Copy(const GHepParticle &particle)
virtual int FirstNonInitStateEntry(void)
virtual GHepParticle * FinalStatePrimaryLepton(void) const
GEvGenMode_t EventGenerationMode(void) const
static int GetPrintLevel()
virtual GHepParticle * TargetNucleus(void) const
virtual vector< int > * GetStableDescendants(int position) const
virtual void ResetRecord(void)
virtual int HitElectronPosition(void) const
virtual GHepParticle * FinalStateHadronicSystem(void) const
bool PolzIsSet(void) const
void GetPolarization(TVector3 &polz)
bool Is2NucleonCluster(int pdgc)
virtual void CompactifyDaughterLists(void)
virtual TBits * EventFlags(void) const
void SetLastMother(int m)
virtual bool IsUnphysical(void) const
double Vz(void) const
Get production z.
ostream & operator<<(ostream &stream, const AlgConfigPool &config_pool)
virtual void Clear(Option_t *opt="")
virtual GHepParticle * HitNucleon(void) const
virtual void SetVertex(double x, double y, double z, double t)
KinePhaseSpace_t fDiffXSecPhSp
specifies which differential cross-section (dsig/dQ2, dsig/dQ2dW, dsig/dxdy,...)
virtual GHepParticle * RemnantNucleus(void) const
virtual void AddParticle(const GHepParticle &p)
double Vy(void) const
Get production y.
double fWeight
event weight
virtual int ProbePosition(void) const
void SetFirstDaughter(int d)
const int kPdgHadronicSyst
void SetUnphysEventMask(const TBits &mask)
GENIE's GHEP MC event record.
Most commonly used PDG codes. A set of utility functions to handle PDG codes is provided in PDGUtils...
STDHEP-like event record entry that can fit a particle or a nucleus.
virtual int FinalStateHadronicSystemPosition(void) const
double fDiffXSec
differential cross section for selected event kinematics
virtual int FinalStatePrimaryLeptonPosition(void) const
virtual int TargetNucleusPosition(void) const
bool IsElectron(int pdgc)
enum genie::EGHepStatus GHepStatus_t
virtual TBits * EventMask(void) const
double Vx(void) const
Get production x.
double Py(void) const
Get Py.