Public Member Functions | Private Member Functions | List of all members
genie::InitialStateAppender Class Reference

Appends the initial state information to the event record. Is a concerete implementation of the EventRecordVisitorI interface. More...

#include <InitialStateAppender.h>

Inheritance diagram for genie::InitialStateAppender:
genie::EventRecordVisitorI genie::Algorithm

Public Member Functions

 InitialStateAppender ()
 
 InitialStateAppender (string config)
 
 ~InitialStateAppender ()
 
void ProcessEventRecord (GHepRecord *event_rec) const
 
- Public Member Functions inherited from genie::EventRecordVisitorI
virtual ~EventRecordVisitorI ()
 
- Public Member Functions inherited from genie::Algorithm
virtual ~Algorithm ()
 
virtual void Configure (const Registry &config)
 
virtual void Configure (string config)
 
virtual void FindConfig (void)
 
virtual const RegistryGetConfig (void) const
 
RegistryGetOwnedConfig (void)
 
virtual const AlgIdId (void) const
 Get algorithm ID. More...
 
virtual AlgStatus_t GetStatus (void) const
 Get algorithm status. More...
 
virtual bool AllowReconfig (void) const
 
virtual AlgCmp_t Compare (const Algorithm *alg) const
 Compare with input algorithm. More...
 
virtual void SetId (const AlgId &id)
 Set algorithm ID. More...
 
virtual void SetId (string name, string config)
 
const AlgorithmSubAlg (const RgKey &registry_key) const
 
void AdoptConfig (void)
 
void AdoptSubstructure (void)
 
virtual void Print (ostream &stream) const
 Print algorithm info. More...
 

Private Member Functions

void AddNeutrino (GHepRecord *event_rec) const
 
void AddNucleus (GHepRecord *event_rec) const
 
void AddStruckParticle (GHepRecord *event_rec) const
 

Additional Inherited Members

- Static Public Member Functions inherited from genie::Algorithm
static string BuildParamVectKey (const std::string &comm_name, unsigned int i)
 
static string BuildParamVectSizeKey (const std::string &comm_name)
 
- Protected Member Functions inherited from genie::EventRecordVisitorI
 EventRecordVisitorI ()
 
 EventRecordVisitorI (string name)
 
 EventRecordVisitorI (string name, string config)
 
- Protected Member Functions inherited from genie::Algorithm
 Algorithm ()
 
 Algorithm (string name)
 
 Algorithm (string name, string config)
 
void Initialize (void)
 
void DeleteConfig (void)
 
void DeleteSubstructure (void)
 
RegistryExtractLocalConfig (const Registry &in) const
 
RegistryExtractLowerConfig (const Registry &in, const string &alg_key) const
 Split an incoming configuration Registry into a block valid for the sub-algo identified by alg_key. More...
 
template<class T >
bool GetParam (const RgKey &name, T &p, bool is_top_call=true) const
 
template<class T >
bool GetParamDef (const RgKey &name, T &p, const T &def) const
 
template<class T >
int GetParamVect (const std::string &comm_name, std::vector< T > &v, bool is_top_call=true) const
 Handle to load vectors of parameters. More...
 
int GetParamVectKeys (const std::string &comm_name, std::vector< RgKey > &k, bool is_top_call=true) const
 
int AddTopRegistry (Registry *rp, bool owns=true)
 add registry with top priority, also update ownership More...
 
int AddLowRegistry (Registry *rp, bool owns=true)
 add registry with lowest priority, also update ownership More...
 
int MergeTopRegistry (const Registry &r)
 
int AddTopRegisties (const vector< Registry * > &rs, bool owns=false)
 Add registries with top priority, also udated Ownerships. More...
 
- Protected Attributes inherited from genie::Algorithm
bool fAllowReconfig
 
bool fOwnsSubstruc
 true if it owns its substructure (sub-algs,...) More...
 
AlgId fID
 algorithm name and configuration set More...
 
vector< Registry * > fConfVect
 
vector< boolfOwnerships
 ownership for every registry in fConfVect More...
 
AlgStatus_t fStatus
 algorithm execution status More...
 
AlgMapfOwnedSubAlgMp
 local pool for owned sub-algs (taken out of the factory pool) More...
 

Detailed Description

Appends the initial state information to the event record. Is a concerete implementation of the EventRecordVisitorI interface.

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

October 04, 2004

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

Definition at line 26 of file InitialStateAppender.h.

Constructor & Destructor Documentation

InitialStateAppender::InitialStateAppender ( )

Definition at line 26 of file InitialStateAppender.cxx.

26  :
27 EventRecordVisitorI("genie::InitialStateAppender")
28 {
29 
30 }
InitialStateAppender::InitialStateAppender ( string  config)

Definition at line 32 of file InitialStateAppender.cxx.

32  :
33 EventRecordVisitorI("genie::InitialStateAppender", config)
34 {
35 
36 }
static Config * config
Definition: config.cpp:1054
InitialStateAppender::~InitialStateAppender ( )

Definition at line 38 of file InitialStateAppender.cxx.

39 {
40 
41 }

Member Function Documentation

void InitialStateAppender::AddNeutrino ( GHepRecord event_rec) const
private

Definition at line 66 of file InitialStateAppender.cxx.

67 {
68  Interaction * interaction = evrec->Summary();
69  const InitialState & init_state = interaction->InitState();
70 
71  TLorentzVector * p4 = init_state.GetProbeP4(kRfLab);
72  const TLorentzVector v4(0.,0.,0.,0.);
73 
74  int pdgc = init_state.ProbePdg();
75 
76  LOG("ISApp", pINFO) << "Adding neutrino [pdgc = " << pdgc << "]";
77 
78  evrec->AddParticle(pdgc,kIStInitialState, -1,-1,-1,-1, *p4, v4);
79 
80  delete p4;
81 }
Summary information for an interaction.
Definition: Interaction.h:56
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
int ProbePdg(void) const
Definition: InitialState.h:64
#define pINFO
Definition: Messenger.h:62
const InitialState & InitState(void) const
Definition: Interaction.h:69
TLorentzVector * GetProbeP4(RefFrame_t rf=kRfHitNucRest) const
Initial State information.
Definition: InitialState.h:48
void InitialStateAppender::AddNucleus ( GHepRecord event_rec) const
private

Definition at line 83 of file InitialStateAppender.cxx.

84 {
85  Interaction * interaction = evrec->Summary();
86  const InitialState & init_state = interaction->InitState();
87  const ProcessInfo & proc_info = interaction->ProcInfo();
88 
89  bool is_nucleus = init_state.Tgt().IsNucleus();
90  if(!is_nucleus && !proc_info.IsGlashowResonance()) {
91  LOG("ISApp", pINFO)
92  << "Not an interaction with a nuclear target - no nucleus to add";
93  return;
94  }
95  int A = init_state.Tgt().A();
96  int Z = init_state.Tgt().Z();
97  int pdgc = pdg::IonPdgCode(A, Z);
98  double M = PDGLibrary::Instance()->Find(pdgc)->Mass();
99 
100  LOG("ISApp", pINFO)
101  << "Adding nucleus [A = " << A << ", Z = " << Z
102  << ", pdg = " << pdgc << "]";
103 
104  evrec->AddParticle(pdgc,kIStInitialState,-1,-1,-1,-1, 0,0,0,M, 0,0,0,0);
105 }
int A(void) const
Definition: Target.h:70
bool IsNucleus(void) const
Definition: Target.cxx:272
Summary information for an interaction.
Definition: Interaction.h:56
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
A class encapsulating an enumeration of interaction types (EM, Weak-CC, Weak-NC) and scattering types...
Definition: ProcessInfo.h:46
int Z(void) const
Definition: Target.h:68
#define pINFO
Definition: Messenger.h:62
static PDGLibrary * Instance(void)
Definition: PDGLibrary.cxx:57
int IonPdgCode(int A, int Z)
Definition: PDGUtils.cxx:68
#define A
Definition: memgrp.cpp:38
const InitialState & InitState(void) const
Definition: Interaction.h:69
const ProcessInfo & ProcInfo(void) const
Definition: Interaction.h:70
TParticlePDG * Find(int pdgc, bool must_exist=true)
Definition: PDGLibrary.cxx:75
const Target & Tgt(void) const
Definition: InitialState.h:66
bool IsGlashowResonance(void) const
Initial State information.
Definition: InitialState.h:48
void InitialStateAppender::AddStruckParticle ( GHepRecord event_rec) const
private

Definition at line 107 of file InitialStateAppender.cxx.

108 {
109  Interaction * interaction = evrec->Summary();
110  const InitialState & init_state = interaction->InitState();
111  const ProcessInfo & proc_info = interaction->ProcInfo();
112 
113  // EDIT: Add the dark matter scattering off electron here
114  bool hit_e = proc_info.IsInverseMuDecay() ||
115  proc_info.IsIMDAnnihilation() ||
116  proc_info.IsNuElectronElastic() ||
117  proc_info.IsDarkMatterElectronElastic() ||
118  proc_info.IsGlashowResonance();
119 
120  if(hit_e) {
121  int pdgc = kPdgElectron;
122  double mass = PDGLibrary::Instance()->Find(pdgc)->Mass();
123  const TLorentzVector p4(0,0,0, mass);
124  const TLorentzVector v4(0.,0.,0.,0.);
125 
126  LOG("ISApp", pINFO) << "Adding struck electron";
127  evrec->AddParticle(pdgc, kIStInitialState, 1, -1, -1, -1, p4, v4);
128  return;
129  }
130 
131  int pdgc = init_state.Tgt().HitNucPdg();
132 
133  if(pdgc != 0) {
134 
135  bool is_nucleus = init_state.Tgt().IsNucleus();
136 
137  GHepStatus_t ist = (is_nucleus) ? kIStNucleonTarget : kIStInitialState;
138  int imom1 = (is_nucleus) ? 1 : -1;
139  int imom2 = -1;
140 
141  const TLorentzVector p4(init_state.Tgt().HitNucP4());
142  const TLorentzVector v4(0.,0.,0.,0.);
143 
144  LOG("ISApp", pINFO)<< "Adding struck nucleon [pdgc = " << pdgc << "]";
145 
146  evrec->AddParticle(pdgc, ist, imom1, imom2, -1, -1, p4, v4);
147 
148  }//if struck nucleon was set
149 }
bool IsDarkMatterElectronElastic(void) const
int HitNucPdg(void) const
Definition: Target.cxx:304
bool IsInverseMuDecay(void) const
bool IsNucleus(void) const
Definition: Target.cxx:272
bool IsIMDAnnihilation(void) const
const int kPdgElectron
Definition: PDGCodes.h:35
Summary information for an interaction.
Definition: Interaction.h:56
const TLorentzVector & HitNucP4(void) const
Definition: Target.h:91
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
bool IsNuElectronElastic(void) const
A class encapsulating an enumeration of interaction types (EM, Weak-CC, Weak-NC) and scattering types...
Definition: ProcessInfo.h:46
#define pINFO
Definition: Messenger.h:62
static PDGLibrary * Instance(void)
Definition: PDGLibrary.cxx:57
const InitialState & InitState(void) const
Definition: Interaction.h:69
const ProcessInfo & ProcInfo(void) const
Definition: Interaction.h:70
TParticlePDG * Find(int pdgc, bool must_exist=true)
Definition: PDGLibrary.cxx:75
const Target & Tgt(void) const
Definition: InitialState.h:66
bool IsGlashowResonance(void) const
enum genie::EGHepStatus GHepStatus_t
Initial State information.
Definition: InitialState.h:48
void InitialStateAppender::ProcessEventRecord ( GHepRecord event_rec) const
virtual

Implements genie::EventRecordVisitorI.

Definition at line 43 of file InitialStateAppender.cxx.

44 {
45 // Adds the initial state particles at the event record (the order is
46 // significant)
47 
48  LOG("ISApp", pINFO) << "Adding the initial state to the event record";
49 
50  //-- add the incoming neutrino to the event record
51  this->AddNeutrino(evrec);
52 
53  //-- add the nuclear target at the event record (if any)
54  this->AddNucleus(evrec);
55 
56  //-- add the struck nucleon to the event record (if any)
57  // It is added with status-code = 0 (init state) if the target was a
58  // free nucleon, or with a status-code = 11 (nucleon target) if the
59  // target was a nucleus.
60  // If the interaction was ve- elastic, inverse muon decay or Glashow
61  // resonance then it will add the target e- instead.
62  this->AddStruckParticle(evrec);
63  LOG("ISApp", pINFO) << *evrec;
64 }
void AddNeutrino(GHepRecord *event_rec) const
void AddNucleus(GHepRecord *event_rec) const
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
void AddStruckParticle(GHepRecord *event_rec) const
#define pINFO
Definition: Messenger.h:62

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