GENIEGenerator
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Public Member Functions | Private Member Functions | Private Attributes | List of all members
genie::DISHadronicSystemGenerator Class Reference

Generates the final state hadronic system in v DIS interactions. Is a concrete implementation of the EventRecordVisitorI interface. More...

#include <DISHadronicSystemGenerator.h>

Inheritance diagram for genie::DISHadronicSystemGenerator:
Inheritance graph
[legend]
Collaboration diagram for genie::DISHadronicSystemGenerator:
Collaboration graph
[legend]

Public Member Functions

 DISHadronicSystemGenerator ()
 
 DISHadronicSystemGenerator (string config)
 
 ~DISHadronicSystemGenerator ()
 
void ProcessEventRecord (GHepRecord *event_rec) const
 
void Configure (const Registry &config)
 
void Configure (string config)
 
- Public Member Functions inherited from genie::HadronicSystemGenerator
void AddTargetNucleusRemnant (GHepRecord *event_rec) const
 
void AddFinalHadronicSyst (GHepRecord *event_rec) const
 
void PreHadronTransportDecays (GHepRecord *event_rec) const
 
TLorentzVector Hadronic4pLAB (GHepRecord *event_rec) const
 
TLorentzVector MomentumTransferLAB (GHepRecord *event_rec) const
 
TVector3 HCM2LAB (GHepRecord *event_rec) const
 
int HadronShowerCharge (GHepRecord *event_rec) const
 
int ResonanceCharge (GHepRecord *event_rec) const
 
- Public Member Functions inherited from genie::EventRecordVisitorI
virtual ~EventRecordVisitorI ()
 
- Public Member Functions inherited from genie::Algorithm
virtual ~Algorithm ()
 
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 SimulateFormationZone (GHepRecord *event_rec) const
 
void LoadConfig (void)
 

Private Attributes

const EventRecordVisitorIfHadronizationModel
 
bool fFilterPreFragmEntries
 
double fR0
 param controling nuclear size More...
 
double fNR
 how far beyond the nuclear boundary does the particle tracker goes? More...
 
double fct0pion
 formation zone (c * formation time) - for pions More...
 
double fct0nucleon
 formation zone (c * formation time) - for nucleons More...
 
double fK
 param multiplying pT^2 in formation zone calculation More...
 

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)
 
static string BuildParamMatKey (const std::string &comm_name, unsigned int i, unsigned int j)
 
static string BuildParamMatRowSizeKey (const std::string &comm_name)
 
static string BuildParamMatColSizeKey (const std::string &comm_name)
 
- Protected Member Functions inherited from genie::HadronicSystemGenerator
 HadronicSystemGenerator ()
 
 HadronicSystemGenerator (string name)
 
 HadronicSystemGenerator (string name, string config)
 
 ~HadronicSystemGenerator ()
 
- 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
 
template<class T >
int GetParamMat (const std::string &comm_name, TMatrixT< T > &mat, bool is_top_call=true) const
 Handle to load matrix of parameters. More...
 
template<class T >
int GetParamMatSym (const std::string &comm_name, TMatrixTSym< T > &mat, bool is_top_call=true) const
 
int GetParamMatKeys (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::HadronicSystemGenerator
const EventRecordVisitorIfPreINukeDecayer
 
- 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< bool > fOwnerships
 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

Generates the final state hadronic system in v DIS interactions. Is a concrete implementation of the EventRecordVisitorI interface.

Author
Costas Andreopoulos <c.andreopoulos cern.ch> University of Liverpool
Created:
October 03, 2004
License:
Copyright (c) 2003-2024, The GENIE Collaboration For the full text of the license visit http://copyright.genie-mc.org

Definition at line 26 of file DISHadronicSystemGenerator.h.

Constructor & Destructor Documentation

DISHadronicSystemGenerator::DISHadronicSystemGenerator ( )

Definition at line 36 of file DISHadronicSystemGenerator.cxx.

36  :
37 HadronicSystemGenerator("genie::DISHadronicSystemGenerator")
38 {
39 
40 }
DISHadronicSystemGenerator::DISHadronicSystemGenerator ( string  config)

Definition at line 42 of file DISHadronicSystemGenerator.cxx.

42  :
43 HadronicSystemGenerator("genie::DISHadronicSystemGenerator", config)
44 {
45 
46 }
DISHadronicSystemGenerator::~DISHadronicSystemGenerator ( )

Definition at line 48 of file DISHadronicSystemGenerator.cxx.

49 {
50 
51 }

Member Function Documentation

void DISHadronicSystemGenerator::Configure ( const Registry config)
virtual

Configure the algorithm with an external registry The registry is merged with the top level registry if it is owned, Otherwise a copy of it is added with the highest priority

Reimplemented from genie::Algorithm.

Definition at line 161 of file DISHadronicSystemGenerator.cxx.

References genie::Algorithm::Configure(), and LoadConfig().

162 {
163  Algorithm::Configure(config);
164  this->LoadConfig();
165 }
virtual void Configure(const Registry &config)
Definition: Algorithm.cxx:62
void DISHadronicSystemGenerator::Configure ( string  config)
virtual

Configure the algorithm from the AlgoConfigPool based on param_set string given in input An algorithm contains a vector of registries coming from different xml configuration files, which are loaded according a very precise prioriy This methods will load a number registries in order of priority: 1) "Tunable" parameter set from CommonParametes. This is loaded with the highest prioriry and it is designed to be used for tuning procedure Usage not expected from the user. 2) For every string defined in "CommonParame" the corresponding parameter set will be loaded from CommonParameter.xml 3) parameter set specified by the config string and defined in the xml file of the algorithm 4) if config is not "Default" also the Default parameter set from the same xml file will be loaded Effectively this avoids the repetion of a parameter when it is not changed in the requested configuration

Reimplemented from genie::Algorithm.

Definition at line 167 of file DISHadronicSystemGenerator.cxx.

References genie::Algorithm::Configure(), and LoadConfig().

168 {
169  Algorithm::Configure(config);
170  this->LoadConfig();
171 }
virtual void Configure(const Registry &config)
Definition: Algorithm.cxx:62
void DISHadronicSystemGenerator::LoadConfig ( void  )
private

Definition at line 173 of file DISHadronicSystemGenerator.cxx.

References fct0nucleon, fct0pion, fFilterPreFragmEntries, fHadronizationModel, fK, fNR, genie::HadronicSystemGenerator::fPreINukeDecayer, fR0, genie::Algorithm::GetParam(), genie::Algorithm::GetParamDef(), LOG, pDEBUG, and genie::Algorithm::SubAlg().

Referenced by Configure().

174 {
176  fPreINukeDecayer = 0;
177 
178  //-- Get the requested hadronization model
180  dynamic_cast<const EventRecordVisitorI *> (this->SubAlg("Hadronizer"));
181  assert(fHadronizationModel);
182 
183  //-- Handle pre-intranuke decays
185  dynamic_cast<const EventRecordVisitorI *> (this->SubAlg("PreTransportDecayer"));
186  assert(fPreINukeDecayer);
187 
188  //-- Get flag to determine whether we copy all fragmentation record entries
189  // into the GHEP record or just the ones marked with kf=1
190  GetParamDef( "FilterPreFragm", fFilterPreFragmEntries, false ) ;
191 
192  //-- Get parameters controlling the nuclear sizes
193  GetParam( "NUCL-R0", fR0 ) ;
194  GetParam( "NUCL-NR", fNR ) ;
195 
196  //-- Get parameters controlling the formation zone simulation
197  GetParam( "FZONE-ct0pion", fct0pion ) ;
198  GetParam( "FZONE-ct0nucleon",fct0nucleon ) ;
199  GetParam( "FZONE-KPt2", fK ) ;
200 
201  LOG("DISHadronicVtx", pDEBUG) << "ct0pion = " << fct0pion << " fermi";
202  LOG("DISHadronicVtx", pDEBUG) << "ct0nucleon = " << fct0nucleon << " fermi";
203  LOG("DISHadronicVtx", pDEBUG) << "K(pt^2) = " << fK;
204 }
double fct0pion
formation zone (c * formation time) - for pions
Defines the EventRecordVisitorI interface. Concrete implementations of this interface use the &#39;Visito...
const EventRecordVisitorI * fPreINukeDecayer
double fK
param multiplying pT^2 in formation zone calculation
const EventRecordVisitorI * fHadronizationModel
double fNR
how far beyond the nuclear boundary does the particle tracker goes?
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
double fct0nucleon
formation zone (c * formation time) - for nucleons
bool GetParamDef(const RgKey &name, T &p, const T &def) const
bool GetParam(const RgKey &name, T &p, bool is_top_call=true) const
double fR0
param controling nuclear size
#define pDEBUG
Definition: Messenger.h:63
const Algorithm * SubAlg(const RgKey &registry_key) const
Definition: Algorithm.cxx:345
void DISHadronicSystemGenerator::ProcessEventRecord ( GHepRecord event_rec) const
virtual

Implements genie::EventRecordVisitorI.

Definition at line 53 of file DISHadronicSystemGenerator.cxx.

References genie::HadronicSystemGenerator::AddFinalHadronicSyst(), fHadronizationModel, genie::HadronicSystemGenerator::Hadronic4pLAB(), genie::Interaction::KinePtr(), genie::Kinematics::SetW(), SimulateFormationZone(), and genie::GHepRecord::Summary().

54 {
55 // This method generates the final state hadronic system
56 
57  //-- Add an entry for the DIS Pre-Fragm. Hadronic State
58  this->AddFinalHadronicSyst(evrec);
59 
60  // set W in the Event Summary
61  //-- Compute the hadronic system invariant mass
62  TLorentzVector p4Had = this->Hadronic4pLAB(evrec);
63 
64  Interaction * interaction = evrec->Summary();
65  interaction->KinePtr()->SetW( p4Had.M() );
66 
67  //-- Add the fragmentation products
69 
70 
71  //-- Simulate the formation zone if not taken directly from the
72  // hadronization model
73  this->SimulateFormationZone(evrec);
74 }
TLorentzVector Hadronic4pLAB(GHepRecord *event_rec) const
Kinematics * KinePtr(void) const
Definition: Interaction.h:76
void ProcessEventRecord(GHepRecord *event_rec) const
const EventRecordVisitorI * fHadronizationModel
Summary information for an interaction.
Definition: Interaction.h:56
void SetW(double W, bool selected=false)
Definition: Kinematics.cxx:279
void SimulateFormationZone(GHepRecord *event_rec) const
void AddFinalHadronicSyst(GHepRecord *event_rec) const
void DISHadronicSystemGenerator::SimulateFormationZone ( GHepRecord event_rec) const
private

Definition at line 76 of file DISHadronicSystemGenerator.cxx.

References genie::units::A, genie::GHepParticle::A(), epsilon, fct0nucleon, fct0pion, genie::GHepRecord::FinalStateHadronicSystem(), fK, fNR, genie::utils::phys::FormationZone(), fR0, genie::pdg::IsIon(), genie::pdg::IsNucleon(), genie::kIStHadronInTheNucleus, LOG, genie::units::m, genie::GHepParticle::Mass(), genie::GHepParticle::Name(), genie::GHepParticle::P4(), pDEBUG, genie::GHepParticle::Pdg(), pINFO, genie::HadronicSystemGenerator::PreHadronTransportDecays(), genie::GHepParticle::SetPosition(), genie::GHepParticle::Status(), genie::GHepRecord::TargetNucleus(), and genie::GHepParticle::X4().

Referenced by ProcessEventRecord().

78 {
79  LOG("DISHadronicVtx", pDEBUG)
80  << "Simulating formation zone for the DIS hadronic system";
81 
82  GHepParticle * nucltgt = evrec->TargetNucleus();
83  if (!nucltgt) {
84  LOG("DISHadronicVtx", pDEBUG)
85  << "No nuclear target was found - No need to simulate formation zones";
86  return;
87  }
88  // Compute the nuclear radius & how far away a particle is being tracked by
89  // the intranuclear hadron transport
90  assert(nucltgt && pdg::IsIon(nucltgt->Pdg()));
91  double A = nucltgt->A();
92  double R = fR0 * TMath::Power(A, 1./3.);
93  R *= TMath::Max(fNR,1.); // particle is tracked much further outside the nuclear boundary as the density is non-zero
94 
95  // Decay very short living particles so that we can give the formation
96  // zone to the daughters
97  this->PreHadronTransportDecays(evrec);
98 
99  // Get hadronic system's 3-momentum
100  GHepParticle * hadronic_system = evrec->FinalStateHadronicSystem();
101  TVector3 p3hadr = hadronic_system->P4()->Vect(); // (px,py,pz)
102 
103  // Loop over GHEP and set the formation zone to the right particles
104  // Limit the maximum formation zone so that particles escaping the
105  // nucleus are placed right outside...
106 
107  TObjArrayIter piter(evrec);
108  GHepParticle * p = 0;
109  int icurr = -1;
110 
111  while( (p = (GHepParticle *) piter.Next()) )
112  {
113  icurr++;
114 
115  GHepStatus_t ist = p->Status();
116  bool apply_formation_zone = (ist==kIStHadronInTheNucleus);
117 
118  if(!apply_formation_zone) continue;
119 
120  LOG("DISHadronicVtx", pINFO)
121  << "Applying formation-zone to " << p->Name();
122 
123  double m = p->Mass();
124  int pdgc = p->Pdg();
125  const TLorentzVector & p4 = *(p->P4());
126  double ct0=0.;
127  pdg::IsNucleon(pdgc) ? ct0=fct0nucleon : ct0=fct0pion;
128  double fz = phys::FormationZone(m,p4,p3hadr,ct0,fK);
129 
130  //-- Apply the formation zone step
131 
132  double step = fz;
133  TVector3 dr = p->P4()->Vect().Unit(); // unit vector along its direction
134  // double c = kLightSpeed / (units::fm/units::ns); // c in fm/nsec
135  dr.SetMag(step); // spatial step size
136  // double dt = step/c; // temporal step:
137  double dt = 0;
138  TLorentzVector dx4(dr,dt); // 4-vector step
139  TLorentzVector x4new = *(p->X4()) + dx4; // new position
140 
141  //-- If the formation zone was large enough that the particle is now outside
142  // the nucleus make sure that it is not placed further away from the
143  // (max distance particles tracked by intranuclear cascade) + ~2 fm
144  double epsilon = 2; // fm
145  double r = x4new.Vect().Mag(); // fm
146  double rmax = R+epsilon;
147  if(r > rmax) {
148  LOG("DISHadronicVtx", pINFO)
149  << "Particle was stepped too far away (r = " << r << " fm)";
150  LOG("DISHadronicVtx", pINFO)
151  << "Placing it ~2 fm away from the furthermost position tracked "
152  << "by intranuclear cascades (r' = " << rmax << " fm)";
153  double scale = rmax/r;
154  x4new *= scale;
155  }
156 
157  p->SetPosition(x4new);
158  }
159 }
double fct0pion
formation zone (c * formation time) - for pions
const TLorentzVector * P4(void) const
Definition: GHepParticle.h:78
bool IsNucleon(int pdgc)
Definition: PDGUtils.cxx:346
double fK
param multiplying pT^2 in formation zone calculation
double Mass(void) const
Mass that corresponds to the PDG code.
GHepStatus_t Status(void) const
Definition: GHepParticle.h:64
double FormationZone(double m, const TLorentzVector &p, const TVector3 &p3hadr, double ct0, double K)
Definition: PhysUtils.cxx:18
const double epsilon
int Pdg(void) const
Definition: GHepParticle.h:63
double fNR
how far beyond the nuclear boundary does the particle tracker goes?
string Name(void) const
Name that corresponds to the PDG code.
void SetPosition(const TLorentzVector &v4)
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
static constexpr double A
Definition: Units.h:74
double fct0nucleon
formation zone (c * formation time) - for nucleons
#define pINFO
Definition: Messenger.h:62
bool IsIon(int pdgc)
Definition: PDGUtils.cxx:42
const TLorentzVector * X4(void) const
Definition: GHepParticle.h:79
int A(void) const
double fR0
param controling nuclear size
static constexpr double m
Definition: Units.h:71
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 PreHadronTransportDecays(GHepRecord *event_rec) const

Member Data Documentation

double genie::DISHadronicSystemGenerator::fct0nucleon
private

formation zone (c * formation time) - for nucleons

Definition at line 53 of file DISHadronicSystemGenerator.h.

Referenced by LoadConfig(), and SimulateFormationZone().

double genie::DISHadronicSystemGenerator::fct0pion
private

formation zone (c * formation time) - for pions

Definition at line 52 of file DISHadronicSystemGenerator.h.

Referenced by LoadConfig(), and SimulateFormationZone().

bool genie::DISHadronicSystemGenerator::fFilterPreFragmEntries
private

Definition at line 49 of file DISHadronicSystemGenerator.h.

Referenced by LoadConfig().

const EventRecordVisitorI* genie::DISHadronicSystemGenerator::fHadronizationModel
private

Definition at line 47 of file DISHadronicSystemGenerator.h.

Referenced by LoadConfig(), and ProcessEventRecord().

double genie::DISHadronicSystemGenerator::fK
private

param multiplying pT^2 in formation zone calculation

Definition at line 54 of file DISHadronicSystemGenerator.h.

Referenced by LoadConfig(), and SimulateFormationZone().

double genie::DISHadronicSystemGenerator::fNR
private

how far beyond the nuclear boundary does the particle tracker goes?

Definition at line 51 of file DISHadronicSystemGenerator.h.

Referenced by LoadConfig(), and SimulateFormationZone().

double genie::DISHadronicSystemGenerator::fR0
private

param controling nuclear size

Definition at line 50 of file DISHadronicSystemGenerator.h.

Referenced by LoadConfig(), and SimulateFormationZone().


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