Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | Static Protected Attributes | Friends | List of all members
genie::GHepRecord Class Reference

GENIE's GHEP MC event record. More...

#include <GHepRecord.h>

Inheritance diagram for genie::GHepRecord:
genie::EventRecord

Public Member Functions

 GHepRecord ()
 
 GHepRecord (int size)
 
 GHepRecord (const GHepRecord &record)
 
 GHepRecord (TRootIOCtor *)
 
virtual ~GHepRecord ()
 
virtual InteractionSummary (void) const
 
virtual void AttachSummary (Interaction *interaction)
 
virtual void AddParticle (const GHepParticle &p)
 
virtual void AddParticle (int pdg, GHepStatus_t ist, int mom1, int mom2, int dau1, int dau2, const TLorentzVector &p, const TLorentzVector &v)
 
virtual void AddParticle (int pdg, GHepStatus_t ist, int mom1, int mom2, int dau1, int dau2, double px, double py, double pz, double E, double x, double y, double z, double t)
 
virtual GHepParticleParticle (int position) const
 
virtual GHepParticleFindParticle (int pdg, GHepStatus_t ist, int start) const
 
virtual int ParticlePosition (int pdg, GHepStatus_t i, int start=0) const
 
virtual int ParticlePosition (GHepParticle *particle, int start=0) const
 
virtual vector< int > * GetStableDescendants (int position) const
 
GEvGenMode_t EventGenerationMode (void) const
 
virtual GHepParticleProbe (void) const
 
virtual GHepParticleTargetNucleus (void) const
 
virtual GHepParticleRemnantNucleus (void) const
 
virtual GHepParticleHitNucleon (void) const
 
virtual GHepParticleHitElectron (void) const
 
virtual GHepParticleFinalStatePrimaryLepton (void) const
 
virtual GHepParticleFinalStateHadronicSystem (void) const
 
virtual int ProbePosition (void) const
 
virtual int TargetNucleusPosition (void) const
 
virtual int RemnantNucleusPosition (void) const
 
virtual int HitNucleonPosition (void) const
 
virtual int HitElectronPosition (void) const
 
virtual int FinalStatePrimaryLeptonPosition (void) const
 
virtual int FinalStateHadronicSystemPosition (void) const
 
virtual unsigned int NEntries (int pdg, GHepStatus_t ist, int start=0) const
 
virtual unsigned int NEntries (int pdg, int start=0) const
 
virtual TBits * EventFlags (void) const
 
virtual TBits * EventMask (void) const
 
virtual bool IsUnphysical (void) const
 
virtual bool Accept (void) const
 
virtual double Weight (void) const
 
virtual double Probability (void) const
 
virtual double XSec (void) const
 
virtual double DiffXSec (void) const
 
virtual KinePhaseSpace_t DiffXSecVars (void) const
 
virtual void SetWeight (double wght)
 
virtual void SetProbability (double prob)
 
virtual void SetXSec (double xsec)
 
virtual void SetDiffXSec (double xsec, KinePhaseSpace_t ps)
 
virtual TLorentzVector * Vertex (void) const
 
virtual void SetVertex (double x, double y, double z, double t)
 
virtual void SetVertex (const TLorentzVector &vtx)
 
virtual void Copy (const GHepRecord &record)
 
virtual void Clear (Option_t *opt="")
 
virtual void ResetRecord (void)
 
virtual void CompactifyDaughterLists (void)
 
virtual void RemoveIntermediateParticles (void)
 
void SetUnphysEventMask (const TBits &mask)
 
void Print (ostream &stream) const
 

Static Public Member Functions

static void SetPrintLevel (int print_level)
 
static int GetPrintLevel ()
 

Protected Member Functions

void InitRecord (void)
 
void CleanRecord (void)
 
virtual void UpdateDaughterLists (void)
 
virtual bool HasCompactDaughterList (int pos)
 
virtual void SwapParticles (int i, int j)
 
virtual void FinalizeDaughterLists (void)
 
virtual int FirstNonInitStateEntry (void)
 

Protected Attributes

InteractionfInteraction
 attached summary information More...
 
TLorentzVector * fVtx
 vertex in the detector coordinate system More...
 
TBits * fEventFlags
 event flags indicating various pathologies or an unphysical event More...
 
TBits * fEventMask
 an input bit-field mask allowing one to ignore bits set in fEventFlags More...
 
double fWeight
 event weight More...
 
double fProb
 event probability (for given flux neutrino and density-weighted path-length for target element) More...
 
double fXSec
 cross section for selected event More...
 
double fDiffXSec
 differential cross section for selected event kinematics More...
 
KinePhaseSpace_t fDiffXSecPhSp
 specifies which differential cross-section (dsig/dQ2, dsig/dQ2dW, dsig/dxdy,...) More...
 

Static Protected Attributes

static int fPrintLevel
 

Friends

ostream & operator<< (ostream &stream, const GHepRecord &event)
 

Detailed Description

GENIE's GHEP MC event record.

Author
Costas Andreopoulos <constantinos.andreopoulos cern.ch> University of Liverpool & STFC Rutherford Appleton Laboratory

October 1, 2004

Copyright (c) 2003-2020, The GENIE Collaboration For the full text of the license visit http://copyright.genie-mc.org

Definition at line 45 of file GHepRecord.h.

Constructor & Destructor Documentation

GHepRecord::GHepRecord ( )

Definition at line 53 of file GHepRecord.cxx.

53  :
54 TClonesArray("genie::GHepParticle")
55 {
56  this->InitRecord();
57 }
void InitRecord(void)
Definition: GHepRecord.cxx:829
GHepRecord::GHepRecord ( int  size)

Definition at line 59 of file GHepRecord.cxx.

59  :
60 TClonesArray("genie::GHepParticle", size)
61 {
62  this->InitRecord();
63 }
void InitRecord(void)
Definition: GHepRecord.cxx:829
decltype(auto) constexpr size(T &&obj)
ADL-aware version of std::size.
Definition: StdUtils.h:92
GHepRecord::GHepRecord ( const GHepRecord record)

Definition at line 65 of file GHepRecord.cxx.

65  :
66 TClonesArray("genie::GHepParticle", record.GetEntries())
67 {
68  this->InitRecord();
69  this->Copy(record);
70 }
virtual void Copy(const GHepRecord &record)
Definition: GHepRecord.cxx:899
void InitRecord(void)
Definition: GHepRecord.cxx:829
GHepRecord::GHepRecord ( TRootIOCtor *  )

Definition at line 72 of file GHepRecord.cxx.

72  :
73 TClonesArray("genie::GHepParticle"),
74 fInteraction(0),
75 fVtx(0),
76 fEventFlags(0),
77 fEventMask(0),
78 fWeight(0.),
79 fProb(0.),
80 fXSec(0.),
81 fDiffXSec(0.)
82 {
83 
84 }
double fProb
event probability (for given flux neutrino and density-weighted path-length for target element) ...
Definition: GHepRecord.h:179
TLorentzVector * fVtx
vertex in the detector coordinate system
Definition: GHepRecord.h:171
TBits * fEventMask
an input bit-field mask allowing one to ignore bits set in fEventFlags
Definition: GHepRecord.h:175
double fXSec
cross section for selected event
Definition: GHepRecord.h:180
Interaction * fInteraction
attached summary information
Definition: GHepRecord.h:168
TBits * fEventFlags
event flags indicating various pathologies or an unphysical event
Definition: GHepRecord.h:174
double fWeight
event weight
Definition: GHepRecord.h:178
double fDiffXSec
differential cross section for selected event kinematics
Definition: GHepRecord.h:181
GHepRecord::~GHepRecord ( )
virtual

Definition at line 86 of file GHepRecord.cxx.

87 {
88  this->CleanRecord();
89 }
void CleanRecord(void)
Definition: GHepRecord.cxx:858

Member Function Documentation

bool GHepRecord::Accept ( void  ) const
virtual

Definition at line 939 of file GHepRecord.cxx.

940 {
941  TBits flags = *fEventFlags;
942  TBits mask = *fEventMask;
943  TBits bitwiseand = flags & mask;
944  bool accept = (bitwiseand.CountBits() == 0);
945  return accept;
946 }
TBits * fEventMask
an input bit-field mask allowing one to ignore bits set in fEventFlags
Definition: GHepRecord.h:175
TBits * fEventFlags
event flags indicating various pathologies or an unphysical event
Definition: GHepRecord.h:174
void GHepRecord::AddParticle ( const GHepParticle p)
virtual

Definition at line 491 of file GHepRecord.cxx.

492 {
493 // Provides a simplified method for inserting entries in the TClonesArray
494 
495  unsigned int pos = this->GetEntries();
496 
497 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
498  LOG("GHEP", pINFO)
499  << "Adding particle with pdgc = " << p.Pdg() << " at slot = " << pos;
500 #endif
501  new ((*this)[pos]) GHepParticle(p);
502 
503  // Update the mother's daughter list. If the newly inserted particle broke
504  // compactification, then run CompactifyDaughterLists()
505  this->UpdateDaughterLists();
506 }
virtual void UpdateDaughterLists(void)
Definition: GHepRecord.cxx:548
int Pdg(void) const
Definition: GHepParticle.h:63
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
#define pINFO
Definition: Messenger.h:62
STDHEP-like event record entry that can fit a particle or a nucleus.
Definition: GHepParticle.h:39
void GHepRecord::AddParticle ( int  pdg,
GHepStatus_t  ist,
int  mom1,
int  mom2,
int  dau1,
int  dau2,
const TLorentzVector &  p,
const TLorentzVector &  v 
)
virtual

Definition at line 508 of file GHepRecord.cxx.

511 {
512 // Provides a 'simplified' method for inserting entries in the TClonesArray
513 
514  unsigned int pos = this->GetEntries();
515 
516 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
517  LOG("GHEP", pINFO)
518  << "Adding particle with pdgc = " << pdg << " at slot = " << pos;
519 #endif
520  new ((*this)[pos]) GHepParticle(pdg,status, mom1,mom2,dau1,dau2, p, v);
521 
522  // Update the mother's daughter list. If the newly inserted particle broke
523  // compactification, then run CompactifyDaughterLists()
524  this->UpdateDaughterLists();
525 }
virtual void UpdateDaughterLists(void)
Definition: GHepRecord.cxx:548
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
p
Definition: test.py:223
#define pINFO
Definition: Messenger.h:62
STDHEP-like event record entry that can fit a particle or a nucleus.
Definition: GHepParticle.h:39
void GHepRecord::AddParticle ( int  pdg,
GHepStatus_t  ist,
int  mom1,
int  mom2,
int  dau1,
int  dau2,
double  px,
double  py,
double  pz,
double  E,
double  x,
double  y,
double  z,
double  t 
)
virtual

Definition at line 527 of file GHepRecord.cxx.

531 {
532 // Provides a 'simplified' method for inserting entries in the TClonesArray
533 
534  unsigned int pos = this->GetEntries();
535 
536 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
537  LOG("GHEP", pINFO)
538  << "Adding particle with pdgc = " << pdg << " at slot = " << pos;
539 #endif
540  new ( (*this)[pos] ) GHepParticle (
541  pdg, status, mom1, mom2, dau1, dau2, px, py, pz, E, x, y, z, t);
542 
543  // Update the mother's daughter list. If the newly inserted particle broke
544  // compactification, then run CompactifyDaughterLists()
545  this->UpdateDaughterLists();
546 }
virtual void UpdateDaughterLists(void)
Definition: GHepRecord.cxx:548
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
#define pINFO
Definition: Messenger.h:62
list x
Definition: train.py:276
STDHEP-like event record entry that can fit a particle or a nucleus.
Definition: GHepParticle.h:39
void GHepRecord::AttachSummary ( Interaction interaction)
virtual

Definition at line 99 of file GHepRecord.cxx.

100 {
102 }
Interaction * fInteraction
attached summary information
Definition: GHepRecord.h:168
void GHepRecord::CleanRecord ( void  )
protected

Definition at line 858 of file GHepRecord.cxx.

859 {
860 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
861  LOG("GHEP", pDEBUG) << "Cleaning up GHepRecord";
862 #endif
863  this->Clear("C");
864 }
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
virtual void Clear(Option_t *opt="")
Definition: GHepRecord.cxx:875
#define pDEBUG
Definition: Messenger.h:63
void GHepRecord::Clear ( Option_t *  opt = "")
virtual

Definition at line 875 of file GHepRecord.cxx.

876 {
877  if (fInteraction) delete fInteraction;
878  fInteraction=0;
879 
880  if (fVtx) delete fVtx;
881  fVtx=0;
882 
883  if(fEventFlags) delete fEventFlags;
884  fEventFlags=0;
885 
886  if(fEventMask) delete fEventMask;
887  fEventMask=0;
888 
889  TClonesArray::Clear(opt);
890 
891 // if (fInteraction) delete fInteraction;
892 // delete fVtx;
893 //
894 // delete fEventFlags;
895 //
896 // TClonesArray::Clear(opt);
897 }
opt
Definition: train.py:196
TLorentzVector * fVtx
vertex in the detector coordinate system
Definition: GHepRecord.h:171
TBits * fEventMask
an input bit-field mask allowing one to ignore bits set in fEventFlags
Definition: GHepRecord.h:175
Interaction * fInteraction
attached summary information
Definition: GHepRecord.h:168
TBits * fEventFlags
event flags indicating various pathologies or an unphysical event
Definition: GHepRecord.h:174
void GHepRecord::CompactifyDaughterLists ( void  )
virtual

compact

Definition at line 643 of file GHepRecord.cxx.

644 {
645  int n = this->GetEntries();
646  if(n<1) return;
647 
648  int i = this->Particle(n-1)->FirstMother();
649  if(i<0) return;
650 
651  // for(int i=0; i<n; i++) {
652  bool compact = this->HasCompactDaughterList(i);
653 
654 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
655  LOG("GHEP", pNOTICE)
656  << "Particle's " << i << " daughter list is "
657  << ((compact) ? "compact" : "__not__ compact");
658 #endif
659 
660  if(!compact) {
661  GHepParticle * p = this->Particle(i);
662 
663  int dau1 = p->FirstDaughter();
664  int dau2 = p->LastDaughter();
665  int ndau = dau2-dau1+1;
666  int ndp = dau2+1;
667  if(dau1==-1) {ndau=0;}
668 
669  int curr_pos = n-1;
670  while (curr_pos > ndp) {
671  this->SwapParticles(curr_pos,curr_pos-1);
672  curr_pos--;
673  }
674  if(ndau>0) {
675  this->Particle(i)->SetFirstDaughter(dau1);
676  this->Particle(i)->SetLastDaughter(dau2+1);
677  } else {
678  this->Particle(i)->SetFirstDaughter(-1);
679  this->Particle(i)->SetLastDaughter(-1);
680  }
681  this->FinalizeDaughterLists();
682 
683  } //!compact
684 
685 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
686  LOG("GHEP", pINFO)
687  << "Done ompactifying daughter-list for particle at slot: " << i;
688 #endif
689  // }
690 }
virtual GHepParticle * Particle(int position) const
Definition: GHepRecord.cxx:104
virtual void SwapParticles(int i, int j)
Definition: GHepRecord.cxx:744
int FirstDaughter(void) const
Definition: GHepParticle.h:68
virtual void FinalizeDaughterLists(void)
Definition: GHepRecord.cxx:791
int FirstMother(void) const
Definition: GHepParticle.h:66
int LastDaughter(void) const
Definition: GHepParticle.h:69
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
void SetLastDaughter(int d)
Definition: GHepParticle.h:135
std::void_t< T > n
virtual bool HasCompactDaughterList(int pos)
Definition: GHepRecord.cxx:692
p
Definition: test.py:223
#define pINFO
Definition: Messenger.h:62
#define pNOTICE
Definition: Messenger.h:61
void SetFirstDaughter(int d)
Definition: GHepParticle.h:134
STDHEP-like event record entry that can fit a particle or a nucleus.
Definition: GHepParticle.h:39
void GHepRecord::Copy ( const GHepRecord record)
virtual

Definition at line 899 of file GHepRecord.cxx.

900 {
901  // clean up
902  this->ResetRecord();
903 
904  // copy event record entries
905  unsigned int ientry = 0;
906  GHepParticle * p = 0;
907  TIter ghepiter(&record);
908  while ( (p = (GHepParticle *) ghepiter.Next()) )
909  new ( (*this)[ientry++] ) GHepParticle(*p);
910 
911  // copy summary
912  fInteraction = new Interaction( *record.fInteraction );
913 
914  // copy flags & mask
915  *fEventFlags = *(record.EventFlags());
916  *fEventMask = *(record.EventMask());
917 
918  // copy vtx position
919  TLorentzVector * v = record.Vertex();
920  fVtx->SetXYZT(v->X(),v->Y(),v->Z(),v->T());
921 
922  // copy weights & xsecs
923  fWeight = record.fWeight;
924  fProb = record.fProb;
925  fXSec = record.fXSec;
926  fDiffXSec = record.fDiffXSec;
927  fDiffXSecPhSp = record.fDiffXSecPhSp;
928 }
double fProb
event probability (for given flux neutrino and density-weighted path-length for target element) ...
Definition: GHepRecord.h:179
TLorentzVector * fVtx
vertex in the detector coordinate system
Definition: GHepRecord.h:171
TBits * fEventMask
an input bit-field mask allowing one to ignore bits set in fEventFlags
Definition: GHepRecord.h:175
virtual TLorentzVector * Vertex(void) const
Definition: GHepRecord.h:140
double fXSec
cross section for selected event
Definition: GHepRecord.h:180
Interaction * fInteraction
attached summary information
Definition: GHepRecord.h:168
TBits * fEventFlags
event flags indicating various pathologies or an unphysical event
Definition: GHepRecord.h:174
p
Definition: test.py:223
virtual void ResetRecord(void)
Definition: GHepRecord.cxx:866
virtual TBits * EventFlags(void) const
Definition: GHepRecord.h:117
KinePhaseSpace_t fDiffXSecPhSp
specifies which differential cross-section (dsig/dQ2, dsig/dQ2dW, dsig/dxdy,...)
Definition: GHepRecord.h:182
double fWeight
event weight
Definition: GHepRecord.h:178
STDHEP-like event record entry that can fit a particle or a nucleus.
Definition: GHepParticle.h:39
double fDiffXSec
differential cross section for selected event kinematics
Definition: GHepRecord.h:181
virtual TBits * EventMask(void) const
Definition: GHepRecord.h:118
virtual double genie::GHepRecord::DiffXSec ( void  ) const
inlinevirtual

Definition at line 127 of file GHepRecord.h.

127 { return fDiffXSec; }
double fDiffXSec
differential cross section for selected event kinematics
Definition: GHepRecord.h:181
virtual KinePhaseSpace_t genie::GHepRecord::DiffXSecVars ( void  ) const
inlinevirtual

Definition at line 128 of file GHepRecord.h.

128 { return fDiffXSecPhSp; }
KinePhaseSpace_t fDiffXSecPhSp
specifies which differential cross-section (dsig/dQ2, dsig/dQ2dW, dsig/dxdy,...)
Definition: GHepRecord.h:182
virtual TBits* genie::GHepRecord::EventFlags ( void  ) const
inlinevirtual

Definition at line 117 of file GHepRecord.h.

117 { return fEventFlags; }
TBits * fEventFlags
event flags indicating various pathologies or an unphysical event
Definition: GHepRecord.h:174
GEvGenMode_t GHepRecord::EventGenerationMode ( void  ) const

Definition at line 209 of file GHepRecord.cxx.

210 {
211  GHepParticle * p0 = this->Particle(0);
212  if(!p0) return kGMdUnknown;
213  GHepParticle * p1 = this->Particle(1);
214  if(!p1) return kGMdUnknown;
215 
216  int p0pdg = p0->Pdg();
217  GHepStatus_t p0st = p0->Status();
218  int p1pdg = p1->Pdg();
219  GHepStatus_t p1st = p1->Status();
220 
221  // In lepton+nucleon/nucleus mode, the 1st entry in the event record
222  // is a charged or neutral lepton with status code = kIStInitialState
223  if( pdg::IsLepton(p0pdg) && p0st == kIStInitialState )
224  {
225  return kGMdLeptonNucleus;
226  }
227 
228  // In dark matter mode, the 1st entry in the event record is a dark
229  // matter particle
230  if( (pdg::IsDarkMatter(p0pdg) || pdg::IsAntiDarkMatter(p0pdg)) && p0st == kIStInitialState )
231  {
232  return kGMdDarkMatterNucleus;
233  }
234 
235  // In hadron+nucleon/nucleus mode, the 1st entry in the event record
236  // is a hadron with status code = kIStInitialState and the 2nd entry
237  // is a nucleon or nucleus with status code = kIStInitialState
238  if( pdg::IsHadron(p0pdg) && p0st == kIStInitialState )
239  {
240  if( (pdg::IsIon(p1pdg) || pdg::IsNucleon(p1pdg)) && p1st == kIStInitialState)
241  {
242  return kGMdHadronNucleus;
243  }
244  }
245 
246  // As above, with a photon as a probe
247  if( p0pdg == kPdgGamma && p0st == kIStInitialState )
248  {
249  if( (pdg::IsIon(p1pdg) || pdg::IsNucleon(p1pdg)) && p1st == kIStInitialState)
250  {
251  return kGMdPhotonNucleus;
252  }
253  }
254 
255  // In nucleon decay mode,
256  // - [if the decayed nucleon was a bound one] the 1st entry in the event
257  // record is a nucleus with status code = kIStInitialState and the
258  // 2nd entry is a nucleon with code = kIStDecayedState
259  // - [if the decayed nucleon was a free one] the first entry in the event
260  // record is a nucleon with status code = kIStInitialState and it has a
261  // single daughter which is a nucleon with status code = kIStDecayedState.
262 
263  if( pdg::IsIon(p0pdg) && p0st == kIStInitialState &&
264  pdg::IsNucleon(p1pdg) && p1st == kIStDecayedState)
265  {
266  return kGMdNucleonDecay;
267  }
268  if( pdg::IsNucleon(p0pdg) && p0st == kIStInitialState &&
269  pdg::IsNucleon(p1pdg) && p1st == kIStDecayedState)
270  {
271  return kGMdNucleonDecay;
272  }
273 
274  return kGMdUnknown;
275 }
virtual GHepParticle * Particle(int position) const
Definition: GHepRecord.cxx:104
bool IsNucleon(int pdgc)
Definition: PDGUtils.cxx:343
bool IsDarkMatter(int pdgc)
Definition: PDGUtils.cxx:124
GHepStatus_t Status(void) const
Definition: GHepParticle.h:64
int Pdg(void) const
Definition: GHepParticle.h:63
bool IsHadron(int pdgc)
Definition: PDGUtils.cxx:389
bool IsAntiDarkMatter(int pdgc)
Definition: PDGUtils.cxx:130
const int kPdgGamma
Definition: PDGCodes.h:189
bool IsIon(int pdgc)
Definition: PDGUtils.cxx:39
STDHEP-like event record entry that can fit a particle or a nucleus.
Definition: GHepParticle.h:39
bool IsLepton(int pdgc)
Definition: PDGUtils.cxx:83
enum genie::EGHepStatus GHepStatus_t
virtual TBits* genie::GHepRecord::EventMask ( void  ) const
inlinevirtual

Definition at line 118 of file GHepRecord.h.

118 { return fEventMask; }
TBits * fEventMask
an input bit-field mask allowing one to ignore bits set in fEventFlags
Definition: GHepRecord.h:175
void GHepRecord::FinalizeDaughterLists ( void  )
protectedvirtual

Definition at line 791 of file GHepRecord.cxx.

792 {
793 // Update all daughter-lists based on particle 'first mother' field.
794 // To work correctly, the daughter-lists must have been compactified first.
795 
796  GHepParticle * p1 = 0;
797  TIter iter1(this);
798  int i1=0;
799  while( (p1 = (GHepParticle *)iter1.Next()) ) {
800  int dau1 = -1;
801  int dau2 = -1;
802  GHepParticle * p2 = 0;
803  TIter iter2(this);
804  int i2=0;
805  while( (p2 = (GHepParticle *)iter2.Next()) ) {
806 
807  if(p2->FirstMother() == i1) {
808  dau1 = (dau1<0) ? i2 : TMath::Min(dau1,i2);
809  dau2 = (dau2<0) ? i2 : TMath::Max(dau2,i2);
810  }
811  i2++;
812  }
813  i1++;
814  p1 -> SetFirstDaughter (dau1);
815  p1 -> SetLastDaughter (dau2);
816  }
817 }
int FirstMother(void) const
Definition: GHepParticle.h:66
STDHEP-like event record entry that can fit a particle or a nucleus.
Definition: GHepParticle.h:39
GHepParticle * GHepRecord::FinalStateHadronicSystem ( void  ) const
virtual

Definition at line 335 of file GHepRecord.cxx.

336 {
337 // Returns the GHepParticle representing the sum of the DIS pre-fragm f/s
338 // hadronic system, or 0 if it does not exist.
339 
340  int ipos = this->FinalStateHadronicSystemPosition();
341  if(ipos>-1) return this->Particle(ipos);
342  return 0;
343 }
virtual GHepParticle * Particle(int position) const
Definition: GHepRecord.cxx:104
virtual int FinalStateHadronicSystemPosition(void) const
Definition: GHepRecord.cxx:463
int GHepRecord::FinalStateHadronicSystemPosition ( void  ) const
virtual

Definition at line 463 of file GHepRecord.cxx.

464 {
465  return this->ParticlePosition(
467 }
virtual int ParticlePosition(int pdg, GHepStatus_t i, int start=0) const
Definition: GHepRecord.cxx:137
const int kPdgHadronicSyst
Definition: PDGCodes.h:210
GHepParticle * GHepRecord::FinalStatePrimaryLepton ( void  ) const
virtual

Definition at line 326 of file GHepRecord.cxx.

327 {
328 // Returns the GHepParticle representing the final state primary lepton.
329 
330  int ipos = this->FinalStatePrimaryLeptonPosition();
331  if(ipos>-1) return this->Particle(ipos);
332  return 0;
333 }
virtual GHepParticle * Particle(int position) const
Definition: GHepRecord.cxx:104
virtual int FinalStatePrimaryLeptonPosition(void) const
Definition: GHepRecord.cxx:451
int GHepRecord::FinalStatePrimaryLeptonPosition ( void  ) const
virtual

Definition at line 451 of file GHepRecord.cxx.

452 {
453 // Returns the GHEP position GHepParticle representing the final state
454 // primary lepton.
455 
456  GHepParticle * probe = this->Probe();
457  if(!probe) return -1;
458 
459  int ifsl = probe->FirstDaughter();
460  return ifsl;
461 }
int FirstDaughter(void) const
Definition: GHepParticle.h:68
virtual GHepParticle * Probe(void) const
Definition: GHepRecord.cxx:277
STDHEP-like event record entry that can fit a particle or a nucleus.
Definition: GHepParticle.h:39
GHepParticle * GHepRecord::FindParticle ( int  pdg,
GHepStatus_t  ist,
int  start 
) const
virtual

Definition at line 118 of file GHepRecord.cxx.

120 {
121 // Returns the first GHepParticle with the input pdg-code and status
122 // starting from the specified position of the event record.
123 
124  int nentries = this->GetEntries();
125  for(int i = start; i < nentries; i++) {
126  GHepParticle * p = (GHepParticle *) (*this)[i];
127  if(p->Status() == status && p->Pdg() == pdg) return p;
128  }
129 
130  LOG("GHEP", pINFO)
131  << "No particle found with: (pos >= " << start
132  << ", pdg = " << pdg << ", ist = " << status << ")";
133 
134  return 0;
135 }
GHepStatus_t Status(void) const
Definition: GHepParticle.h:64
int Pdg(void) const
Definition: GHepParticle.h:63
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
p
Definition: test.py:223
#define pINFO
Definition: Messenger.h:62
STDHEP-like event record entry that can fit a particle or a nucleus.
Definition: GHepParticle.h:39
int GHepRecord::FirstNonInitStateEntry ( void  )
protectedvirtual

Definition at line 731 of file GHepRecord.cxx.

732 {
733  GHepParticle * p = 0;
734  TIter iter(this);
735  int pos = 0;
736  while( (p = (GHepParticle *)iter.Next()) ) {
737  int ist = p->Status();
738  if(ist != kIStInitialState && ist != kIStNucleonTarget) return pos;
739  pos++;
740  }
741  return pos;
742 }
GHepStatus_t Status(void) const
Definition: GHepParticle.h:64
p
Definition: test.py:223
STDHEP-like event record entry that can fit a particle or a nucleus.
Definition: GHepParticle.h:39
int GHepRecord::GetPrintLevel ( )
static

Definition at line 952 of file GHepRecord.cxx.

953 {
954  return fPrintLevel;
955 }
static int fPrintLevel
Definition: GHepRecord.h:196
vector< int > * GHepRecord::GetStableDescendants ( int  position) const
virtual

Definition at line 174 of file GHepRecord.cxx.

175 {
176 // Returns a list of all stable descendants of the GHEP entry in the input
177 // slot. The user adopts the output vector.
178 
179  // return null if particle index is out of range
180  int nentries = this->GetEntries();
181  if(position<0 || position>=nentries) return 0;
182 
183  vector<int> * descendants = new vector<int>;
184 
185  // return itself if it is a stable final state particle
186  if(this->Particle(position)->Status() == kIStStableFinalState) {
187  descendants->push_back(position);
188  return descendants;
189  }
190 
191  for(int i = 0; i < nentries; i++) {
192  if(i==position) continue;
193  GHepParticle * p = (GHepParticle *) (*this)[i];
194  if(p->Status() != kIStStableFinalState) continue;
195  bool is_descendant=false;
196  int mom = p->FirstMother();
197  while(mom>-1) {
198  if(mom==position) is_descendant=true;
199  if(is_descendant) {
200  descendants->push_back(i);
201  break;
202  }
203  mom = this->Particle(mom)->FirstMother();
204  }
205  }
206  return descendants;
207 }
virtual GHepParticle * Particle(int position) const
Definition: GHepRecord.cxx:104
GHepStatus_t Status(void) const
Definition: GHepParticle.h:64
int FirstMother(void) const
Definition: GHepParticle.h:66
p
Definition: test.py:223
STDHEP-like event record entry that can fit a particle or a nucleus.
Definition: GHepParticle.h:39
bool GHepRecord::HasCompactDaughterList ( int  pos)
protectedvirtual

Definition at line 692 of file GHepRecord.cxx.

693 {
694 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
695  LOG("GHEP", pDEBUG) << "Examining daughter-list of particle at: " << pos;
696 #endif
697  vector<int> daughters;
698  GHepParticle * p = 0;
699  TIter iter(this);
700  int i=0;
701  while( (p = (GHepParticle *)iter.Next()) ) {
702  if(p->FirstMother() == pos) {
703 
704 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
705  LOG("GHEP", pDEBUG) << "Particle at: " << i << " is a daughter";
706 #endif
707  daughters.push_back(i);
708  }
709  i++;
710  }
711 
712  bool is_compact = true;
713  if(daughters.size()>1) {
714  sort(daughters.begin(), daughters.end());
715  vector<int>::iterator diter = daughters.begin();
716  int prev = *diter;
717  for(; diter != daughters.end(); ++diter) {
718  int curr = *diter;
719  is_compact = is_compact && (TMath::Abs(prev-curr)<=1);
720  prev = curr;
721  }
722  }
723 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
724  LOG("GHEP", pINFO)
725  << "Daughter-list of particle at: " << pos << " is "
726  << (is_compact ? "" : "not") << " compact";
727 #endif
728  return is_compact;
729 }
intermediate_table::iterator iterator
int FirstMother(void) const
Definition: GHepParticle.h:66
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
size_t size
Definition: lodepng.cpp:55
p
Definition: test.py:223
#define pINFO
Definition: Messenger.h:62
STDHEP-like event record entry that can fit a particle or a nucleus.
Definition: GHepParticle.h:39
#define pDEBUG
Definition: Messenger.h:63
GHepParticle * GHepRecord::HitElectron ( void  ) const
virtual

Definition at line 316 of file GHepRecord.cxx.

317 {
318 // Returns the GHepParticle representing the struck electron, or 0 if it does
319 // not exist.
320 
321  int ipos = this->HitElectronPosition();
322  if(ipos>-1) return this->Particle(ipos);
323  return 0;
324 }
virtual GHepParticle * Particle(int position) const
Definition: GHepRecord.cxx:104
virtual int HitElectronPosition(void) const
Definition: GHepRecord.cxx:433
int GHepRecord::HitElectronPosition ( void  ) const
virtual

Definition at line 433 of file GHepRecord.cxx.

434 {
435 // Returns the GHEP position of the GHepParticle representing a hit electron.
436 // Same as above..
437 
438  GHepParticle * nucleus = this->TargetNucleus();
439 
440  int ipos = (nucleus) ? 2 : 1;
441 
442  GHepParticle * p = this->Particle(ipos);
443  if(!p) return -1;
444 
445  bool ise = pdg::IsElectron(p->Pdg());
446  if(ise && p->Status()==kIStInitialState) return ipos;
447 
448  return -1;
449 }
virtual GHepParticle * Particle(int position) const
Definition: GHepRecord.cxx:104
GHepStatus_t Status(void) const
Definition: GHepParticle.h:64
int Pdg(void) const
Definition: GHepParticle.h:63
p
Definition: test.py:223
virtual GHepParticle * TargetNucleus(void) const
Definition: GHepRecord.cxx:286
STDHEP-like event record entry that can fit a particle or a nucleus.
Definition: GHepParticle.h:39
bool IsElectron(int pdgc)
Definition: PDGUtils.cxx:185
GHepParticle * GHepRecord::HitNucleon ( void  ) const
virtual

Definition at line 306 of file GHepRecord.cxx.

307 {
308 // Returns the GHepParticle representing the struck nucleon, or 0 if it does
309 // not exist.
310 
311  int ipos = this->HitNucleonPosition();
312  if(ipos>-1) return this->Particle(ipos);
313  return 0;
314 }
virtual GHepParticle * Particle(int position) const
Definition: GHepRecord.cxx:104
virtual int HitNucleonPosition(void) const
Definition: GHepRecord.cxx:410
int GHepRecord::HitNucleonPosition ( void  ) const
virtual

Definition at line 410 of file GHepRecord.cxx.

411 {
412 // Returns the GHEP position of the GHepParticle representing the hit nucleon.
413 // If a struck nucleon is set it will be at slot 2 (for scattering off nuclear
414 // targets) or at slot 1 (for free nucleon scattering).
415 // If the struck nucleon is not set (eg coherent scattering, ve- scattering)
416 // it returns 0.
417 
418  GHepParticle * nucleus = this->TargetNucleus();
419 
420  int ipos = (nucleus) ? 2 : 1;
421  GHepStatus_t ist = (nucleus) ? kIStNucleonTarget : kIStInitialState;
422 
423  GHepParticle * p = this->Particle(ipos);
424  if(!p) return -1;
425 
426 // bool isN = pdg::IsNeutronOrProton(p->Pdg());
427  bool isN = pdg::IsNucleon(p->Pdg()) || pdg::Is2NucleonCluster(p->Pdg());
428  if(isN && p->Status()==ist) return ipos;
429 
430  return -1;
431 }
virtual GHepParticle * Particle(int position) const
Definition: GHepRecord.cxx:104
bool IsNucleon(int pdgc)
Definition: PDGUtils.cxx:343
GHepStatus_t Status(void) const
Definition: GHepParticle.h:64
int Pdg(void) const
Definition: GHepParticle.h:63
p
Definition: test.py:223
virtual GHepParticle * TargetNucleus(void) const
Definition: GHepRecord.cxx:286
bool Is2NucleonCluster(int pdgc)
Definition: PDGUtils.cxx:399
STDHEP-like event record entry that can fit a particle or a nucleus.
Definition: GHepParticle.h:39
enum genie::EGHepStatus GHepStatus_t
void GHepRecord::InitRecord ( void  )
protected

Definition at line 829 of file GHepRecord.cxx.

830 {
831 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
832  LOG("GHEP", pDEBUG) << "Initializing GHepRecord";
833 #endif
834  fInteraction = 0;
835  fWeight = 1.;
836  fProb = 1.;
837  fXSec = 0.;
838  fDiffXSec = 0.;
840  fVtx = new TLorentzVector(0,0,0,0);
841 
842  fEventFlags = new TBits(GHepFlags::NFlags());
843  fEventFlags -> ResetAllBits(false);
844 
845  fEventMask = new TBits(GHepFlags::NFlags());
846 //fEventMask -> ResetAllBits(true);
847  for(unsigned int i = 0; i < GHepFlags::NFlags(); i++) {
848  fEventMask->SetBitNumber(i, true);
849  }
850 
851  LOG("GHEP", pINFO)
852  << "Initialised unphysical event mask (bits: " << GHepFlags::NFlags() - 1
853  << " -> 0) : " << *fEventMask;
854 
855  this->SetOwner(true);
856 }
double fProb
event probability (for given flux neutrino and density-weighted path-length for target element) ...
Definition: GHepRecord.h:179
TLorentzVector * fVtx
vertex in the detector coordinate system
Definition: GHepRecord.h:171
TBits * fEventMask
an input bit-field mask allowing one to ignore bits set in fEventFlags
Definition: GHepRecord.h:175
double fXSec
cross section for selected event
Definition: GHepRecord.h:180
Interaction * fInteraction
attached summary information
Definition: GHepRecord.h:168
TBits * fEventFlags
event flags indicating various pathologies or an unphysical event
Definition: GHepRecord.h:174
static unsigned int NFlags(void)
Definition: GHepFlags.h:76
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
#define pINFO
Definition: Messenger.h:62
KinePhaseSpace_t fDiffXSecPhSp
specifies which differential cross-section (dsig/dQ2, dsig/dQ2dW, dsig/dxdy,...)
Definition: GHepRecord.h:182
double fWeight
event weight
Definition: GHepRecord.h:178
double fDiffXSec
differential cross section for selected event kinematics
Definition: GHepRecord.h:181
#define pDEBUG
Definition: Messenger.h:63
virtual bool genie::GHepRecord::IsUnphysical ( void  ) const
inlinevirtual

Definition at line 119 of file GHepRecord.h.

119 { return (fEventFlags->CountBits()>0); }
TBits * fEventFlags
event flags indicating various pathologies or an unphysical event
Definition: GHepRecord.h:174
unsigned int GHepRecord::NEntries ( int  pdg,
GHepStatus_t  ist,
int  start = 0 
) const
virtual

Definition at line 469 of file GHepRecord.cxx.

470 {
471  unsigned int nentries = 0;
472 
473  for(int i = start; i < this->GetEntries(); i++) {
474  GHepParticle * p = (GHepParticle *) (*this)[i];
475  if(p->Pdg()==pdg && p->Status()==ist) nentries++;
476  }
477  return nentries;
478 }
GHepStatus_t Status(void) const
Definition: GHepParticle.h:64
int Pdg(void) const
Definition: GHepParticle.h:63
p
Definition: test.py:223
STDHEP-like event record entry that can fit a particle or a nucleus.
Definition: GHepParticle.h:39
unsigned int GHepRecord::NEntries ( int  pdg,
int  start = 0 
) const
virtual

Definition at line 480 of file GHepRecord.cxx.

481 {
482  unsigned int nentries = 0;
483 
484  for(int i = start; i < this->GetEntries(); i++) {
485  GHepParticle * p = (GHepParticle *) (*this)[i];
486  if(p->Pdg()==pdg) nentries++;
487  }
488  return nentries;
489 }
int Pdg(void) const
Definition: GHepParticle.h:63
p
Definition: test.py:223
STDHEP-like event record entry that can fit a particle or a nucleus.
Definition: GHepParticle.h:39
GHepParticle * GHepRecord::Particle ( int  position) const
virtual

Definition at line 104 of file GHepRecord.cxx.

105 {
106 // Returns the GHepParticle from the specified position of the event record.
107 
108  if( position >=0 && position < this->GetEntries() ) {
109  GHepParticle * particle = (GHepParticle *) (*this)[position];
110  if(particle) return particle;
111  }
112  LOG("GHEP", pINFO)
113  << "No particle found with: (pos = " << position << ")";
114 
115  return 0;
116 }
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
#define pINFO
Definition: Messenger.h:62
STDHEP-like event record entry that can fit a particle or a nucleus.
Definition: GHepParticle.h:39
int GHepRecord::ParticlePosition ( int  pdg,
GHepStatus_t  i,
int  start = 0 
) const
virtual

Definition at line 137 of file GHepRecord.cxx.

139 {
140 // Returns the position of the first GHepParticle with the input pdg-code
141 // and status starting from the specified position of the event record.
142 
143  int nentries = this->GetEntries();
144  for(int i = start; i < nentries; i++) {
145  GHepParticle * p = (GHepParticle *) (*this)[i];
146  if(p->Status() == status && p->Pdg() == pdg) return i;
147  }
148 
149  LOG("GHEP", pINFO)
150  << "No particle found with: (pos >= " << start
151  << ", pdg = " << pdg << ", ist = " << status << ")";
152 
153  return -1;
154 }
GHepStatus_t Status(void) const
Definition: GHepParticle.h:64
int Pdg(void) const
Definition: GHepParticle.h:63
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
p
Definition: test.py:223
#define pINFO
Definition: Messenger.h:62
STDHEP-like event record entry that can fit a particle or a nucleus.
Definition: GHepParticle.h:39
int GHepRecord::ParticlePosition ( GHepParticle particle,
int  start = 0 
) const
virtual

Definition at line 156 of file GHepRecord.cxx.

157 {
158 // Returns the position of the first match with the specified GHepParticle
159 // starting from the specified position of the event record.
160 
161  int nentries = this->GetEntries();
162  for(int i = start; i < nentries; i++) {
163  GHepParticle * p = (GHepParticle *) (*this)[i];
164  if( p->Compare(particle) ) return i;
165  }
166 
167  LOG("GHEP", pINFO)
168  << "No particle found with pos >= " << start
169  << " matching particle: " << *particle;
170 
171  return -1;
172 }
bool Compare(const GHepParticle *p) const
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
p
Definition: test.py:223
#define pINFO
Definition: Messenger.h:62
STDHEP-like event record entry that can fit a particle or a nucleus.
Definition: GHepParticle.h:39
void GHepRecord::Print ( ostream &  stream) const

Definition at line 957 of file GHepRecord.cxx.

958 {
959  // Print levels:
960  // 0 -> prints particle list
961  // 1 -> prints particle list + event flags
962  // 2 -> prints particle list + event flags + wght/xsec
963  // 3 -> prints particle list + event flags + wght/xsec + summary
964  // 10 -> as in level 0 but showing particle positions too
965  // 11 -> as in level 1 but showing particle positions too
966  // 12 -> as in level 2 but showing particle positions too
967  // 13 -> as in level 3 but showing particle positions too
968 
969  bool accept_input_print_level =
970  (fPrintLevel >= 0 && fPrintLevel <= 3) ||
971  (fPrintLevel >=10 && fPrintLevel <=13);
972 
973  int printlevel = (accept_input_print_level) ? fPrintLevel : 3;
974  int printlevel_orig = printlevel;
975 
976  bool showpos = false;
977  if(printlevel >= 10) {
978  printlevel-=10;
979  showpos=true;
980  }
981 
982  stream << "\n\n|";
983  stream << setfill('-') << setw(115) << "|";
984 
985  stream << "\n|GENIE GHEP Event Record [print level: "
986  << setfill(' ') << setw(3) << printlevel_orig << "]"
987  << setfill(' ') << setw(73) << "|";
988 
989  stream << "\n|";
990  stream << setfill('-') << setw(115) << "|";
991 
992  stream << "\n| ";
993  stream << setfill(' ') << setw(6) << "Idx | "
994  << setfill(' ') << setw(16) << "Name | "
995  << setfill(' ') << setw(6) << "Ist | "
996  << setfill(' ') << setw(13) << "PDG | "
997  << setfill(' ') << setw(12) << "Mother | "
998  << setfill(' ') << setw(12) << "Daughter | "
999  << setfill(' ') << setw(10) << ((showpos) ? "Px(x) |" : "Px | ")
1000  << setfill(' ') << setw(10) << ((showpos) ? "Py(y) |" : "Py | ")
1001  << setfill(' ') << setw(10) << ((showpos) ? "Pz(z) |" : "Pz | ")
1002  << setfill(' ') << setw(10) << ((showpos) ? "E(t) |" : "E | ")
1003  << setfill(' ') << setw(10) << "m | ";
1004 
1005  stream << "\n|";
1006  stream << setfill('-') << setw(115) << "|";
1007 
1008  GHepParticle * p = 0;
1009  TObjArrayIter piter(this);
1010  TVector3 polarization(0,0,0);
1011 
1012  unsigned int idx = 0;
1013 
1014  double sum_E = 0;
1015  double sum_px = 0;
1016  double sum_py = 0;
1017  double sum_pz = 0;
1018 
1019  while( (p = (GHepParticle *) piter.Next()) ) {
1020 
1021  stream << "\n| ";
1022  stream << setfill(' ') << setw(3) << idx++ << " | ";
1023  stream << setfill(' ') << setw(13) << p->Name() << " | ";
1024  stream << setfill(' ') << setw(3) << p->Status() << " | ";
1025  stream << setfill(' ') << setw(10) << p->Pdg() << " | ";
1026  stream << setfill(' ') << setw(3) << p->FirstMother() << " | ";
1027  stream << setfill(' ') << setw(3) << p->LastMother() << " | ";
1028  stream << setfill(' ') << setw(3) << p->FirstDaughter() << " | ";
1029  stream << setfill(' ') << setw(3) << p->LastDaughter() << " | ";
1030  stream << std::fixed << setprecision(3);
1031  stream << setfill(' ') << setw(7) << p->Px() << " | ";
1032  stream << setfill(' ') << setw(7) << p->Py() << " | ";
1033  stream << setfill(' ') << setw(7) << p->Pz() << " | ";
1034  stream << setfill(' ') << setw(7) << p->E() << " | ";
1035 
1036  if (p->IsOnMassShell())
1037  stream << setfill(' ') << setw(7) << p->Mass() << " | ";
1038  else
1039  stream << setfill('*') << setw(7) << p->Mass() << " | M = "
1040  << p->P4()->M() << " ";
1041 
1042  if (p->PolzIsSet()) {
1043  p->GetPolarization(polarization);
1044  stream << "P = (" << polarization.x() << "," << polarization.y()
1045  << "," << polarization.z() << ")";
1046  }
1047 
1048  if (p->RescatterCode() != -1) {
1049  stream << "FSI = " << p->RescatterCode();
1050  }
1051 
1052  // plot particle position if requested
1053  if(showpos) {
1054  stream << "\n| ";
1055  stream << setfill(' ') << setw(6) << " | ";
1056  stream << setfill(' ') << setw(16) << " | ";
1057  stream << setfill(' ') << setw(6) << " | ";
1058  stream << setfill(' ') << setw(13) << " | ";
1059  stream << setfill(' ') << setw(6) << " | ";
1060  stream << setfill(' ') << setw(6) << " | ";
1061  stream << setfill(' ') << setw(6) << " | ";
1062  stream << setfill(' ') << setw(6) << " | ";
1063  stream << std::fixed << setprecision(3);
1064  stream << setfill(' ') << setw(7) << p->Vx() << " | ";
1065  stream << setfill(' ') << setw(7) << p->Vy() << " | ";
1066  stream << setfill(' ') << setw(7) << p->Vz() << " | ";
1067  stream << setfill(' ') << setw(7) << p->Vt() << " | ";
1068  stream << setfill(' ') << setw(10) << " | ";
1069  }
1070 
1071  // compute P4_{final} - P4_{nitial}
1072 
1073  if(p->Status() == kIStStableFinalState ||
1075 
1076  sum_E += p->E();
1077  sum_px += p->Px();
1078  sum_py += p->Py();
1079  sum_pz += p->Pz();
1080  } else
1081  if(p->Status() == kIStInitialState) {
1082  /*
1083  if(p->Status() == kIStInitialState || p->Status() == kIStNucleonTarget) {
1084  */
1085  sum_E -= p->E();
1086  sum_px -= p->Px();
1087  sum_py -= p->Py();
1088  sum_pz -= p->Pz();
1089  }
1090 
1091  } // loop over particles
1092 
1093  stream << "\n|";
1094  stream << setfill('-') << setw(115) << "|";
1095 
1096  // Print SUMS
1097  stream << "\n| ";
1098  stream << setfill(' ') << setw(17) << "Fin-Init: "
1099  << setfill(' ') << setw(6) << " "
1100  << setfill(' ') << setw(18) << " "
1101  << setfill(' ') << setw(12) << " "
1102  << setfill(' ') << setw(12) << " | ";
1103  stream << std::fixed << setprecision(3);
1104  stream << setfill(' ') << setw(7) << sum_px << " | ";
1105  stream << setfill(' ') << setw(7) << sum_py << " | ";
1106  stream << setfill(' ') << setw(7) << sum_pz << " | ";
1107  stream << setfill(' ') << setw(7) << sum_E << " | ";
1108  stream << setfill(' ') << setw(10) << " | ";
1109 
1110  stream << "\n|";
1111  stream << setfill('-') << setw(115) << "|";
1112 
1113  // Print vertex
1114 
1115  GHepParticle * probe = this->Probe();
1116  if(probe){
1117  stream << "\n| ";
1118  stream << setfill(' ') << setw(17) << "Vertex: ";
1119  stream << setfill(' ') << setw(11)
1120  << ((probe) ? probe->Name() : "unknown probe") << " @ (";
1121 
1122  stream << std::fixed << setprecision(5);
1123  stream << "x = " << setfill(' ') << setw(11) << this->Vertex()->X() << " m, ";
1124  stream << "y = " << setfill(' ') << setw(11) << this->Vertex()->Y() << " m, ";
1125  stream << "z = " << setfill(' ') << setw(11) << this->Vertex()->Z() << " m, ";
1126  stream << std::scientific << setprecision(6);
1127  stream << "t = " << setfill(' ') << setw(15) << this->Vertex()->T() << " s) ";
1128  stream << std::fixed << setprecision(3);
1129  stream << setfill(' ') << setw(2) << "|";
1130 
1131  stream << "\n|";
1132  stream << setfill('-') << setw(115) << "|";
1133  }
1134 
1135  // Print FLAGS
1136 
1137  if(printlevel>=1) {
1138  stream << "\n| ";
1139  stream << "Err flag [bits:" << fEventFlags->GetNbits()-1 << "->0] : "
1140  << *fEventFlags << " | "
1141  << "1st set: " << setfill(' ') << setw(56)
1142  << ( this->IsUnphysical() ?
1143  GHepFlags::Describe(GHepFlag_t(fEventFlags->FirstSetBit())) :
1144  "none") << " | ";
1145  stream << "\n| ";
1146  stream << "Err mask [bits:" << fEventMask->GetNbits()-1 << "->0] : "
1147  << *fEventMask << " | "
1148  << "Is unphysical: " << setfill(' ') << setw(5)
1149  << utils::print::BoolAsYNString(this->IsUnphysical()) << " | "
1150  << "Accepted: " << setfill(' ') << setw(5)
1152  << " |";
1153  stream << "\n|";
1154  stream << setfill('-') << setw(115) << "|";
1155  }
1156 
1157  if(printlevel>=2) {
1158  stream << "\n| ";
1159  stream << std::scientific << setprecision(5);
1160 
1161  stream << "sig(Ev) = "
1162  << setfill(' ') << setw(17) << fXSec/units::cm2
1163  << " cm^2 |";
1164 
1165  switch(fDiffXSecPhSp) {
1166  case ( kPSyfE ) :
1167  stream << " dsig(y;E)/dy = " << setfill(' ') << setw(13) << fDiffXSec/units::cm2 << " cm^2 |";
1168  break;
1169  case ( kPSxyfE ) :
1170  stream << " d2sig(x,y;E)/dxdy = " << setfill(' ') << setw(13) << fDiffXSec/units::cm2 << " cm^2 |";
1171  break;
1172  case ( kPSxytfE ) :
1173  stream << " d3sig(x,y,t;E)/dxdydt = " << setfill(' ') << setw(13) << fDiffXSec/units::cm2 << " cm^2/GeV^2 |";
1174  break;
1175  case ( kPSQ2fE ) :
1176  stream << " dsig(Q2;E)/dQ2 = " << setfill(' ') << setw(13) << fDiffXSec/units::cm2 << " cm^2/GeV^2 |";
1177  break;
1178  case ( kPSQ2vfE ) :
1179  stream << " dsig(Q2,v;E)/dQ2dv = " << setfill(' ') << setw(13) << fDiffXSec/units::cm2 << " cm^2/GeV^3 |";
1180  break;
1181  case ( kPSWQ2fE ) :
1182  stream << " d2sig(W,Q2;E)/dWdQ2 = " << setfill(' ') << setw(13) << fDiffXSec/units::cm2 << " cm^2/GeV^3 |";
1183  break;
1184  default :
1185  stream << " dsig(Ev;{K_s})/dK = " << setfill(' ') << setw(13) << fDiffXSec/units::cm2 << " cm^2/{K} |";
1186  }
1187  stream << " Weight = "
1188  << setfill(' ') << setw(16)
1189  << std::fixed << setprecision(5)
1190  << fWeight
1191  << " |";
1192 
1193  stream << "\n|";
1194  stream << setfill('-') << setw(115) << "|";
1195  }
1196 
1197  stream << "\n";
1198  stream << setfill(' ');
1199 
1200  if(printlevel==3) {
1202  else stream << "NULL Interaction!" << endl;
1203  }
1204  stream << "\n";
1205 }
int RescatterCode(void) const
Definition: GHepParticle.h:65
bool IsOnMassShell(void) const
double E(void) const
Get energy.
Definition: GHepParticle.h:91
const TLorentzVector * P4(void) const
Definition: GHepParticle.h:78
int FirstDaughter(void) const
Definition: GHepParticle.h:68
string BoolAsYNString(bool b)
Definition: PrintUtils.cxx:108
enum genie::EGHepFlag GHepFlag_t
TBits * fEventMask
an input bit-field mask allowing one to ignore bits set in fEventFlags
Definition: GHepRecord.h:175
double Mass(void) const
Mass that corresponds to the PDG code.
double Pz(void) const
Get Pz.
Definition: GHepParticle.h:90
static const char * Describe(GHepFlag_t flag)
Definition: GHepFlags.h:42
GHepStatus_t Status(void) const
Definition: GHepParticle.h:64
Q_EXPORT QTSManip setprecision(int p)
Definition: qtextstream.h:343
virtual TLorentzVector * Vertex(void) const
Definition: GHepRecord.h:140
double Px(void) const
Get Px.
Definition: GHepParticle.h:88
virtual GHepParticle * Probe(void) const
Definition: GHepRecord.cxx:277
int LastMother(void) const
Definition: GHepParticle.h:67
double Vt(void) const
Get production time.
Definition: GHepParticle.h:97
int Pdg(void) const
Definition: GHepParticle.h:63
int FirstMother(void) const
Definition: GHepParticle.h:66
string Name(void) const
Name that corresponds to the PDG code.
double fXSec
cross section for selected event
Definition: GHepRecord.h:180
Interaction * fInteraction
attached summary information
Definition: GHepRecord.h:168
int LastDaughter(void) const
Definition: GHepParticle.h:69
TBits * fEventFlags
event flags indicating various pathologies or an unphysical event
Definition: GHepRecord.h:174
virtual bool Accept(void) const
Definition: GHepRecord.cxx:939
static constexpr double cm2
Definition: Units.h:69
p
Definition: test.py:223
bool PolzIsSet(void) const
Q_EXPORT QTSManip setw(int w)
Definition: qtextstream.h:331
void GetPolarization(TVector3 &polz)
virtual bool IsUnphysical(void) const
Definition: GHepRecord.h:119
double Vz(void) const
Get production z.
Definition: GHepParticle.h:96
KinePhaseSpace_t fDiffXSecPhSp
specifies which differential cross-section (dsig/dQ2, dsig/dQ2dW, dsig/dxdy,...)
Definition: GHepRecord.h:182
double Vy(void) const
Get production y.
Definition: GHepParticle.h:95
double fWeight
event weight
Definition: GHepRecord.h:178
Q_EXPORT QTSManip setfill(int f)
Definition: qtextstream.h:337
STDHEP-like event record entry that can fit a particle or a nucleus.
Definition: GHepParticle.h:39
double fDiffXSec
differential cross section for selected event kinematics
Definition: GHepRecord.h:181
static int fPrintLevel
Definition: GHepRecord.h:196
QTextStream & endl(QTextStream &s)
double Vx(void) const
Get production x.
Definition: GHepParticle.h:94
double Py(void) const
Get Py.
Definition: GHepParticle.h:89
virtual double genie::GHepRecord::Probability ( void  ) const
inlinevirtual

Definition at line 125 of file GHepRecord.h.

125 { return fProb; }
double fProb
event probability (for given flux neutrino and density-weighted path-length for target element) ...
Definition: GHepRecord.h:179
GHepParticle * GHepRecord::Probe ( void  ) const
virtual

Definition at line 277 of file GHepRecord.cxx.

278 {
279 // Returns the GHepParticle representing the probe (neutrino, e,...).
280 
281  int ipos = this->ProbePosition();
282  if(ipos>-1) return this->Particle(ipos);
283  return 0;
284 }
virtual GHepParticle * Particle(int position) const
Definition: GHepRecord.cxx:104
virtual int ProbePosition(void) const
Definition: GHepRecord.cxx:345
int GHepRecord::ProbePosition ( void  ) const
virtual

Definition at line 345 of file GHepRecord.cxx.

346 {
347 // Returns the GHEP position of the GHepParticle representing the probe
348 // (neutrino, e,...).
349 
350  // The probe is *always* at slot 0.
351  GEvGenMode_t mode = this->EventGenerationMode();
352  if(mode == kGMdLeptonNucleus ||
353  mode == kGMdDarkMatterNucleus ||
354  mode == kGMdHadronNucleus ||
355  mode == kGMdPhotonNucleus)
356  {
357  return 0;
358  }
359  return -1;
360 }
Enumeration of GENIE event generation modes.
GEvGenMode_t EventGenerationMode(void) const
Definition: GHepRecord.cxx:209
GHepParticle * GHepRecord::RemnantNucleus ( void  ) const
virtual

Definition at line 296 of file GHepRecord.cxx.

297 {
298 // Returns the GHepParticle representing the remnant nucleus,
299 // or 0 if it does not exist.
300 
301  int ipos = this->RemnantNucleusPosition();
302  if(ipos>-1) return this->Particle(ipos);
303  return 0;
304 }
virtual GHepParticle * Particle(int position) const
Definition: GHepRecord.cxx:104
virtual int RemnantNucleusPosition(void) const
Definition: GHepRecord.cxx:389
int GHepRecord::RemnantNucleusPosition ( void  ) const
virtual

Definition at line 389 of file GHepRecord.cxx.

390 {
391 // Returns the GHEP position of the GHepParticle representing the remnant
392 // nucleus - or -1 if the interaction takes place at a free nucleon.
393 
394  GHepParticle * p = this->TargetNucleus();
395  if(!p) return -1;
396 
397  int dau1 = p->FirstDaughter();
398  int dau2 = p->LastDaughter();
399 
400  if(dau1==-1 && dau2==-1) return -1;
401 
402  for(int i=dau1; i<=dau2; i++) {
403  GHepParticle * dp = this->Particle(i);
404  int dpdgc = dp->Pdg();
405  if(pdg::IsIon(dpdgc) && dp->Status()==kIStStableFinalState) return i;
406  }
407  return -1;
408 }
virtual GHepParticle * Particle(int position) const
Definition: GHepRecord.cxx:104
int FirstDaughter(void) const
Definition: GHepParticle.h:68
GHepStatus_t Status(void) const
Definition: GHepParticle.h:64
int Pdg(void) const
Definition: GHepParticle.h:63
int LastDaughter(void) const
Definition: GHepParticle.h:69
p
Definition: test.py:223
virtual GHepParticle * TargetNucleus(void) const
Definition: GHepRecord.cxx:286
bool IsIon(int pdgc)
Definition: PDGUtils.cxx:39
STDHEP-like event record entry that can fit a particle or a nucleus.
Definition: GHepParticle.h:39
void GHepRecord::RemoveIntermediateParticles ( void  )
virtual

Definition at line 609 of file GHepRecord.cxx.

610 {
611  LOG("GHEP", pNOTICE) << "Removing all intermediate particles from GHEP";
612  this->Compress();
613 
614  int i=0;
615  GHepParticle * p = 0;
616 
617  TIter iter(this);
618  while( (p = (GHepParticle *)iter.Next()) ) {
619 
620  if(!p) continue;
621  GHepStatus_t ist = p->Status();
622 
623  bool keep = (ist==kIStInitialState) ||
624  (ist==kIStStableFinalState) || (ist==kIStNucleonTarget);
625  if(keep) {
626  p->SetFirstDaughter(-1);
627  p->SetLastDaughter(-1);
628  p->SetFirstMother(-1);
629  p->SetLastMother(-1);
630  } else {
631  LOG("GHEP", pNOTICE)
632  << "Removing: " << p->Name() << " from slot: " << i;
633  this->RemoveAt(i);
634  }
635  i++;
636  }
637 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
638  LOG("GHEP", pDEBUG) << "Compressing GHEP record to remove empty slots";
639 #endif
640  this->Compress();
641 }
void SetFirstMother(int m)
Definition: GHepParticle.h:132
GHepStatus_t Status(void) const
Definition: GHepParticle.h:64
string Name(void) const
Name that corresponds to the PDG code.
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
void SetLastDaughter(int d)
Definition: GHepParticle.h:135
p
Definition: test.py:223
void SetLastMother(int m)
Definition: GHepParticle.h:133
#define pNOTICE
Definition: Messenger.h:61
void SetFirstDaughter(int d)
Definition: GHepParticle.h:134
void Compress(gar::raw::ADCvector_t &adc, gar::raw::Compress_t compress)
In-place compression of raw data buffer.
Definition: raw.cxx:23
STDHEP-like event record entry that can fit a particle or a nucleus.
Definition: GHepParticle.h:39
enum genie::EGHepStatus GHepStatus_t
#define pDEBUG
Definition: Messenger.h:63
void GHepRecord::ResetRecord ( void  )
virtual

Definition at line 866 of file GHepRecord.cxx.

867 {
868 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
869  LOG("GHEP", pDEBUG) << "Reseting GHepRecord";
870 #endif
871  this->CleanRecord();
872  this->InitRecord();
873 }
void InitRecord(void)
Definition: GHepRecord.cxx:829
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
#define pDEBUG
Definition: Messenger.h:63
void CleanRecord(void)
Definition: GHepRecord.cxx:858
virtual void genie::GHepRecord::SetDiffXSec ( double  xsec,
KinePhaseSpace_t  ps 
)
inlinevirtual

Definition at line 133 of file GHepRecord.h.

134  { fDiffXSecPhSp = ps;
135  fDiffXSec = (xsec>0) ? xsec : 0.;
136  }
static constexpr double ps
Definition: Units.h:99
KinePhaseSpace_t fDiffXSecPhSp
specifies which differential cross-section (dsig/dQ2, dsig/dQ2dW, dsig/dxdy,...)
Definition: GHepRecord.h:182
double fDiffXSec
differential cross section for selected event kinematics
Definition: GHepRecord.h:181
void GHepRecord::SetPrintLevel ( int  print_level)
static

Definition at line 948 of file GHepRecord.cxx.

949 {
950  fPrintLevel = print_level;
951 }
static int fPrintLevel
Definition: GHepRecord.h:196
virtual void genie::GHepRecord::SetProbability ( double  prob)
inlinevirtual

Definition at line 131 of file GHepRecord.h.

131 { fProb = (prob>0) ? prob : 0.; }
double fProb
event probability (for given flux neutrino and density-weighted path-length for target element) ...
Definition: GHepRecord.h:179
void GHepRecord::SetUnphysEventMask ( const TBits &  mask)

Definition at line 930 of file GHepRecord.cxx.

931 {
932  *fEventMask = mask;
933 
934  LOG("GHEP", pINFO)
935  << "Setting unphysical event mask (bits: " << GHepFlags::NFlags() - 1
936  << " -> 0) : " << *fEventMask;
937 }
TBits * fEventMask
an input bit-field mask allowing one to ignore bits set in fEventFlags
Definition: GHepRecord.h:175
static unsigned int NFlags(void)
Definition: GHepFlags.h:76
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
#define pINFO
Definition: Messenger.h:62
void GHepRecord::SetVertex ( double  x,
double  y,
double  z,
double  t 
)
virtual

Definition at line 819 of file GHepRecord.cxx.

820 {
821  fVtx->SetXYZT(x,y,z,t);
822 }
TLorentzVector * fVtx
vertex in the detector coordinate system
Definition: GHepRecord.h:171
list x
Definition: train.py:276
void GHepRecord::SetVertex ( const TLorentzVector &  vtx)
virtual

Definition at line 824 of file GHepRecord.cxx.

825 {
826  fVtx->SetXYZT(vtx.X(),vtx.Y(),vtx.Z(),vtx.T());
827 }
TLorentzVector * fVtx
vertex in the detector coordinate system
Definition: GHepRecord.h:171
virtual void genie::GHepRecord::SetWeight ( double  wght)
inlinevirtual

Definition at line 130 of file GHepRecord.h.

130 { fWeight = (wght>0) ? wght : 0.; }
double fWeight
event weight
Definition: GHepRecord.h:178
virtual void genie::GHepRecord::SetXSec ( double  xsec)
inlinevirtual

Definition at line 132 of file GHepRecord.h.

132 { fXSec = (xsec>0) ? xsec : 0.; }
double fXSec
cross section for selected event
Definition: GHepRecord.h:180
Interaction * GHepRecord::Summary ( void  ) const
virtual

Definition at line 91 of file GHepRecord.cxx.

92 {
93  if(!fInteraction) {
94  LOG("GHEP", pWARN) << "Returning NULL interaction";
95  }
96  return fInteraction;
97 }
Interaction * fInteraction
attached summary information
Definition: GHepRecord.h:168
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
#define pWARN
Definition: Messenger.h:60
void GHepRecord::SwapParticles ( int  i,
int  j 
)
protectedvirtual

Definition at line 744 of file GHepRecord.cxx.

745 {
746 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
747  LOG("GHEP", pINFO) << "Swapping GHepParticles : " << i << " <--> " << j;
748 #endif
749  int n = this->GetEntries();
750  assert(i>=0 && j>=0 && i<n && j<n);
751 
752  if(i==j) return;
753 
754  GHepParticle * pi = this->Particle(i);
755  GHepParticle * pj = this->Particle(j);
756  GHepParticle * tmp = new GHepParticle(*pi);
757 
758  pi->Copy(*pj);
759  pj->Copy(*tmp);
760 
761  delete tmp;
762 
763  // tell their daughters
764  if(pi->HasDaughters()) {
765 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
766  LOG("GHEP", pINFO)
767  << pi->Name() << "(previously at pos: " << j
768  << ") is now at pos: " << i << " -> Notify daughters";
769 #endif
770  for(int k=0; k<n; k++) {
771  if(this->Particle(k)->FirstMother()==j) {
772  this->Particle(k)->SetFirstMother(i);
773  }
774  }
775  }
776 
777  if(pj->HasDaughters()) {
778 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
779  LOG("GHEP", pINFO)
780  << pj->Name() << "(previously at pos: " << i
781  << ") is now at pos: " << j << " -> Notify daughters";
782 #endif
783  for(int k=0; k<n; k++) {
784  if(this->Particle(k)->FirstMother()==i) {
785  this->Particle(k)->SetFirstMother(j);
786  }
787  }
788  }
789 }
void SetFirstMother(int m)
Definition: GHepParticle.h:132
virtual GHepParticle * Particle(int position) const
Definition: GHepRecord.cxx:104
string Name(void) const
Name that corresponds to the PDG code.
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
bool HasDaughters(void) const
Definition: GHepParticle.h:70
std::void_t< T > n
void Copy(const GHepParticle &particle)
#define pINFO
Definition: Messenger.h:62
string tmp
Definition: languages.py:63
float pi
Definition: units.py:11
STDHEP-like event record entry that can fit a particle or a nucleus.
Definition: GHepParticle.h:39
GHepParticle * GHepRecord::TargetNucleus ( void  ) const
virtual

Definition at line 286 of file GHepRecord.cxx.

287 {
288 // Returns the GHepParticle representing the target / initial state nucleus,
289 // or 0 if it does not exist.
290 
291  int ipos = this->TargetNucleusPosition();
292  if(ipos>-1) return this->Particle(ipos);
293  return 0;
294 }
virtual GHepParticle * Particle(int position) const
Definition: GHepRecord.cxx:104
virtual int TargetNucleusPosition(void) const
Definition: GHepRecord.cxx:362
int GHepRecord::TargetNucleusPosition ( void  ) const
virtual

Definition at line 362 of file GHepRecord.cxx.

363 {
364 // Returns the GHEP position of the GHepParticle representing the target
365 // nucleus - or -1 if the interaction takes place at a free nucleon.
366 
367  GEvGenMode_t mode = this->EventGenerationMode();
368 
369  if(mode == kGMdLeptonNucleus ||
370  mode == kGMdDarkMatterNucleus ||
371  mode == kGMdHadronNucleus ||
372  mode == kGMdPhotonNucleus)
373  {
374  GHepParticle * p = this->Particle(1); // If exists, it will be at slot 1
375  if(!p) return -1;
376  int pdgc = p->Pdg();
377  if(pdg::IsIon(pdgc) && p->Status()==kIStInitialState) return 1;
378  }
379  if(mode == kGMdNucleonDecay) {
380  GHepParticle * p = this->Particle(0); // If exists, it will be at slot 0
381  if(!p) return -1;
382  int pdgc = p->Pdg();
383  if(pdg::IsIon(pdgc) && p->Status()==kIStInitialState) return 0;
384  }
385 
386  return -1;
387 }
virtual GHepParticle * Particle(int position) const
Definition: GHepRecord.cxx:104
GHepStatus_t Status(void) const
Definition: GHepParticle.h:64
int Pdg(void) const
Definition: GHepParticle.h:63
Enumeration of GENIE event generation modes.
GEvGenMode_t EventGenerationMode(void) const
Definition: GHepRecord.cxx:209
p
Definition: test.py:223
bool IsIon(int pdgc)
Definition: PDGUtils.cxx:39
STDHEP-like event record entry that can fit a particle or a nucleus.
Definition: GHepParticle.h:39
void GHepRecord::UpdateDaughterLists ( void  )
protectedvirtual

Definition at line 548 of file GHepRecord.cxx.

549 {
550  int pos = this->GetEntries() - 1; // position of last entry
551 
552 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
553  LOG("GHEP", pINFO)
554  << "Updating the daughter-list for the mother of particle at: " << pos;
555 #endif
556 
557  GHepParticle * p = this->Particle(pos);
558  assert(p);
559 
560  int mom_pos = p->FirstMother();
561 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
562  LOG("GHEP", pINFO) << "Mother particle is at slot: " << mom_pos;
563 #endif
564  if(mom_pos==-1) return; // may not have mom (eg init state)
565  GHepParticle * mom = this->Particle(mom_pos);
566  if(!mom) return; // may not have mom (eg init state)
567 
568  int dau1 = mom->FirstDaughter();
569  int dau2 = mom->LastDaughter();
570 
571  // handles the case where the daughter list was initially empty
572  if(dau1 == -1) {
573  mom->SetFirstDaughter(pos);
574  mom->SetLastDaughter(pos);
575  #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
576  LOG("GHEP", pINFO)
577  << "Done! Daughter-list is compact: [" << pos << ", " << pos << "]";
578  #endif
579  return;
580  }
581  // handles the case where the new daughter is added at the slot just before
582  // an already compact daughter list
583  if(pos == dau1-1) {
584  mom->SetFirstDaughter(pos);
585  #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
586  LOG("GHEP", pINFO)
587  << "Done! Daughter-list is compact: [" << pos << ", " << dau2 << "]";
588  #endif
589  return;
590  }
591  // handles the case where the new daughter is added at the slot just after
592  // an already compact daughter list
593  if(pos == dau2+1) {
594  mom->SetLastDaughter(pos);
595  #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
596  LOG("GHEP", pINFO)
597  << "Done! Daughter-list is compact: [" << dau1 << ", " << pos << "]";
598  #endif
599  return;
600  }
601 
602  // If you are here, then the last particle insertion broke the daughter
603  // list compactification - Run the compactifier
604  LOG("GHEP", pNOTICE)
605  << "Daughter-list is not compact - Running compactifier";
606  this->CompactifyDaughterLists();
607 }
virtual GHepParticle * Particle(int position) const
Definition: GHepRecord.cxx:104
int FirstDaughter(void) const
Definition: GHepParticle.h:68
int FirstMother(void) const
Definition: GHepParticle.h:66
int LastDaughter(void) const
Definition: GHepParticle.h:69
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
void SetLastDaughter(int d)
Definition: GHepParticle.h:135
p
Definition: test.py:223
#define pINFO
Definition: Messenger.h:62
virtual void CompactifyDaughterLists(void)
Definition: GHepRecord.cxx:643
#define pNOTICE
Definition: Messenger.h:61
void SetFirstDaughter(int d)
Definition: GHepParticle.h:134
STDHEP-like event record entry that can fit a particle or a nucleus.
Definition: GHepParticle.h:39
virtual TLorentzVector* genie::GHepRecord::Vertex ( void  ) const
inlinevirtual

Definition at line 140 of file GHepRecord.h.

140 { return fVtx; }
TLorentzVector * fVtx
vertex in the detector coordinate system
Definition: GHepRecord.h:171
virtual double genie::GHepRecord::Weight ( void  ) const
inlinevirtual

Definition at line 124 of file GHepRecord.h.

124 { return fWeight; }
double fWeight
event weight
Definition: GHepRecord.h:178
virtual double genie::GHepRecord::XSec ( void  ) const
inlinevirtual

Definition at line 126 of file GHepRecord.h.

126 { return fXSec; }
double fXSec
cross section for selected event
Definition: GHepRecord.h:180

Friends And Related Function Documentation

ostream& operator<< ( ostream &  stream,
const GHepRecord event 
)
friend

Definition at line 46 of file GHepRecord.cxx.

47  {
48  rec.Print(stream);
49  return stream;
50  }
rec
Definition: tracks.py:88

Member Data Documentation

double genie::GHepRecord::fDiffXSec
protected

differential cross section for selected event kinematics

Definition at line 181 of file GHepRecord.h.

KinePhaseSpace_t genie::GHepRecord::fDiffXSecPhSp
protected

specifies which differential cross-section (dsig/dQ2, dsig/dQ2dW, dsig/dxdy,...)

Definition at line 182 of file GHepRecord.h.

TBits* genie::GHepRecord::fEventFlags
protected

event flags indicating various pathologies or an unphysical event

Definition at line 174 of file GHepRecord.h.

TBits* genie::GHepRecord::fEventMask
protected

an input bit-field mask allowing one to ignore bits set in fEventFlags

Definition at line 175 of file GHepRecord.h.

Interaction* genie::GHepRecord::fInteraction
protected

attached summary information

Definition at line 168 of file GHepRecord.h.

int genie::GHepRecord::fPrintLevel
staticprotected

Definition at line 196 of file GHepRecord.h.

double genie::GHepRecord::fProb
protected

event probability (for given flux neutrino and density-weighted path-length for target element)

Definition at line 179 of file GHepRecord.h.

TLorentzVector* genie::GHepRecord::fVtx
protected

vertex in the detector coordinate system

Definition at line 171 of file GHepRecord.h.

double genie::GHepRecord::fWeight
protected

event weight

Definition at line 178 of file GHepRecord.h.

double genie::GHepRecord::fXSec
protected

cross section for selected event

Definition at line 180 of file GHepRecord.h.


The documentation for this class was generated from the following files: