Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4TheoFSGenerator Class Reference

#include <G4TheoFSGenerator.hh>

Inheritance diagram for G4TheoFSGenerator:
G4HadronicInteraction

Public Member Functions

 G4TheoFSGenerator (const G4String &name="TheoFSGenerator")
 
 ~G4TheoFSGenerator ()
 
G4HadFinalStateApplyYourself (const G4HadProjectile &thePrimary, G4Nucleus &theNucleus)
 
void SetTransport (G4VIntraNuclearTransportModel *const value)
 
void SetHighEnergyGenerator (G4VHighEnergyGenerator *const value)
 
void SetQuasiElasticChannel (G4QuasiElasticChannel *const value)
 
virtual std::pair< G4double,
G4double
GetEnergyMomentumCheckLevels () const
 
void ModelDescription (std::ostream &outFile) const
 
- Public Member Functions inherited from G4HadronicInteraction
 G4HadronicInteraction (const G4String &modelName="HadronicModel")
 
virtual ~G4HadronicInteraction ()
 
virtual G4double SampleInvariantT (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
 
virtual G4bool IsApplicable (const G4HadProjectile &, G4Nucleus &)
 
G4double GetMinEnergy () const
 
G4double GetMinEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMinEnergy (G4double anEnergy)
 
void SetMinEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMinEnergy (G4double anEnergy, const G4Material *aMaterial)
 
G4double GetMaxEnergy () const
 
G4double GetMaxEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMaxEnergy (const G4double anEnergy)
 
void SetMaxEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMaxEnergy (G4double anEnergy, const G4Material *aMaterial)
 
const G4HadronicInteractionGetMyPointer () const
 
virtual G4int GetVerboseLevel () const
 
virtual void SetVerboseLevel (G4int value)
 
const G4StringGetModelName () const
 
void DeActivateFor (const G4Material *aMaterial)
 
void ActivateFor (const G4Material *aMaterial)
 
void DeActivateFor (const G4Element *anElement)
 
void ActivateFor (const G4Element *anElement)
 
G4bool IsBlocked (const G4Material *aMaterial) const
 
G4bool IsBlocked (const G4Element *anElement) const
 
void SetRecoilEnergyThreshold (G4double val)
 
G4double GetRecoilEnergyThreshold () const
 
G4bool operator== (const G4HadronicInteraction &right) const
 
G4bool operator!= (const G4HadronicInteraction &right) const
 
virtual const std::pair
< G4double, G4double
GetFatalEnergyCheckLevels () const
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 

Additional Inherited Members

- Protected Member Functions inherited from G4HadronicInteraction
void SetModelName (const G4String &nam)
 
G4bool IsBlocked () const
 
void Block ()
 
- Protected Attributes inherited from G4HadronicInteraction
G4HadFinalState theParticleChange
 
G4int verboseLevel
 
G4double theMinEnergy
 
G4double theMaxEnergy
 
G4bool isBlocked
 

Detailed Description

Definition at line 53 of file G4TheoFSGenerator.hh.

Constructor & Destructor Documentation

G4TheoFSGenerator::G4TheoFSGenerator ( const G4String name = "TheoFSGenerator")

Definition at line 40 of file G4TheoFSGenerator.cc.

41  : G4HadronicInteraction(name)
42  , theTransport(0), theHighEnergyGenerator(0)
43  , theQuasielastic(0)
44  {
45  theParticleChange = new G4HadFinalState;
46 }
G4HadronicInteraction(const G4String &modelName="HadronicModel")
G4TheoFSGenerator::~G4TheoFSGenerator ( )

Definition at line 48 of file G4TheoFSGenerator.cc.

49 {
50  delete theParticleChange;
51 }

Member Function Documentation

G4HadFinalState * G4TheoFSGenerator::ApplyYourself ( const G4HadProjectile thePrimary,
G4Nucleus theNucleus 
)
virtual

Implements G4HadronicInteraction.

Definition at line 68 of file G4TheoFSGenerator.cc.

References G4HadFinalState::AddSecondary(), G4HadFinalState::Clear(), G4cout, G4endl, G4UniformRand, G4HadProjectile::Get4Momentum(), G4Nucleus::GetA_asInt(), G4HadProjectile::GetDefinition(), G4QuasiElasticChannel::GetFraction(), G4IonTable::GetIonMass(), G4ParticleTable::GetIonTable(), G4HadProjectile::GetKineticEnergy(), G4V3DNucleus::GetMassNumber(), G4V3DNucleus::GetNucleons(), G4ParticleTable::GetParticleTable(), G4ParticleDefinition::GetPDGMass(), G4VHighEnergyGenerator::GetProjectileNucleus(), G4VHighEnergyGenerator::GetWoundedNucleus(), G4Nucleus::GetZ_asInt(), isAlive, G4Neutron::Neutron(), G4DecayStrongResonances::Propagate(), G4VIntraNuclearTransportModel::Propagate(), G4VIntraNuclearTransportModel::PropagateNuclNucl(), G4Proton::Proton(), G4QuasiElasticChannel::Scatter(), G4VHighEnergyGenerator::Scatter(), G4HadFinalState::SetEnergyChange(), G4HadFinalState::SetMomentumChange(), G4VIntraNuclearTransportModel::SetPrimaryProjectile(), G4HadFinalState::SetStatusChange(), stopAndKill, CLHEP::Hep3Vector::unit(), and CLHEP::HepLorentzVector::vect().

69 {
70  // init particle change
71  theParticleChange->Clear();
72  theParticleChange->SetStatusChange(stopAndKill);
73 
74  // check if models have been registered, and use default, in case this is not true @@
75 
76  // get result from high energy model
77  G4DynamicParticle aTemp(const_cast<G4ParticleDefinition *>(thePrimary.GetDefinition()),
78  thePrimary.Get4Momentum().vect());
79  const G4DynamicParticle * aPart = &aTemp;
80 
81  if ( theQuasielastic ) {
82 
83  if ( theQuasielastic->GetFraction(theNucleus, *aPart) > G4UniformRand() )
84  {
85  //G4cout<<"___G4TheoFSGenerator: before Scatter (1) QE=" << theQuasielastic<<G4endl;
86  G4KineticTrackVector *result= theQuasielastic->Scatter(theNucleus, *aPart);
87  //G4cout << "^^G4TheoFSGenerator: after Scatter (1) " << G4endl;
88  if (result)
89  {
90  for(unsigned int i=0; i<result->size(); i++)
91  {
92  G4DynamicParticle * aNew =
93  new G4DynamicParticle(result->operator[](i)->GetDefinition(),
94  result->operator[](i)->Get4Momentum().e(),
95  result->operator[](i)->Get4Momentum().vect());
96  theParticleChange->AddSecondary(aNew);
97  delete result->operator[](i);
98  }
99  delete result;
100 
101  } else
102  {
103  theParticleChange->SetStatusChange(isAlive);
104  theParticleChange->SetEnergyChange(thePrimary.GetKineticEnergy());
105  theParticleChange->SetMomentumChange(thePrimary.Get4Momentum().vect().unit());
106 
107  }
108  return theParticleChange;
109  }
110  }
111 
112  G4KineticTrackVector * theInitialResult =
113  theHighEnergyGenerator->Scatter(theNucleus, *aPart);
114 
115 //#define DEBUG_initial_result
116  #ifdef DEBUG_initial_result
117  G4double E_out(0);
119  std::vector<G4KineticTrack *>::iterator ir_iter;
120  for(ir_iter=theInitialResult->begin(); ir_iter!=theInitialResult->end(); ir_iter++)
121  {
122  //G4cout << "TheoFS secondary, mom " << (*ir_iter)->GetDefinition()->GetParticleName() << " " << (*ir_iter)->Get4Momentum() << G4endl;
123  E_out += (*ir_iter)->Get4Momentum().e();
124  }
125  G4double init_mass= ionTable->GetIonMass(theNucleus.GetZ_asInt(),theNucleus.GetA_asInt());
126  G4double init_E=aPart->Get4Momentum().e();
127  // residual nucleus
128 
129  const std::vector<G4Nucleon> & thy = theHighEnergyGenerator->GetWoundedNucleus()->GetNucleons();
130 
131  G4int resZ(0),resA(0);
132  G4double delta_m(0);
133  for(size_t them=0; them<thy.size(); them++)
134  {
135  if(thy[them].AreYouHit()) {
136  ++resA;
137  if ( thy[them].GetDefinition() == G4Proton::Proton() ) {
138  ++resZ;
139  delta_m +=G4Proton::Proton()->GetPDGMass();
140  } else {
141  delta_m +=G4Neutron::Neutron()->GetPDGMass();
142  }
143  }
144  }
145 
146  G4double final_mass(0);
147  if ( theNucleus.GetA_asInt() ) {
148  final_mass=ionTable->GetIonMass(theNucleus.GetZ_asInt()-resZ,theNucleus.GetA_asInt()- resA);
149  }
150  G4double E_excit=init_mass + init_E - final_mass - E_out;
151  G4cout << " Corrected delta mass " << init_mass - final_mass - delta_m << G4endl;
152  G4cout << "initial E, mass = " << init_E << ", " << init_mass << G4endl;
153  G4cout << " final E, mass = " << E_out <<", " << final_mass << " excitation_E " << E_excit << G4endl;
154  #endif
155 
156  G4ReactionProductVector * theTransportResult = NULL;
157 
158 // Uzhi Nov. 2012
159  G4V3DNucleus* theProjectileNucleus = theHighEnergyGenerator->GetProjectileNucleus();
160 if(theProjectileNucleus == 0) // Uzhi Nov. 2012
161 { // Uzhi Nov. 2012
162 
163  G4int hitCount = 0;
164  const std::vector<G4Nucleon>& they = theHighEnergyGenerator->GetWoundedNucleus()->GetNucleons();
165  for(size_t them=0; them<they.size(); them++)
166  {
167  if(they[them].AreYouHit()) hitCount ++;
168  }
169  if(hitCount != theHighEnergyGenerator->GetWoundedNucleus()->GetMassNumber() )
170  {
171  theTransport->SetPrimaryProjectile(thePrimary); // For Bertini Cascade
172  theTransportResult =
173  theTransport->Propagate(theInitialResult, theHighEnergyGenerator->GetWoundedNucleus());
174  if ( !theTransportResult ) {
175  G4cout << "G4TheoFSGenerator: null ptr from transport propagate " << G4endl;
176  throw G4HadronicException(__FILE__, __LINE__, "Null ptr from transport propagate");
177  }
178  }
179  else
180  {
181  theTransportResult = theDecay.Propagate(theInitialResult, theHighEnergyGenerator->GetWoundedNucleus());
182  if ( !theTransportResult ) {
183  G4cout << "G4TheoFSGenerator: null ptr from decay propagate " << G4endl;
184  throw G4HadronicException(__FILE__, __LINE__, "Null ptr from decay propagate");
185  }
186  }
187 
188 } else // Uzhi Nov. 2012
189 { // Uzhi Nov. 2012
190  theTransport->SetPrimaryProjectile(thePrimary);
191  theTransportResult =
192  theTransport->PropagateNuclNucl(theInitialResult,
193  theHighEnergyGenerator->GetWoundedNucleus(),
194  theHighEnergyGenerator->GetProjectileNucleus());
195  if ( !theTransportResult ) {
196  G4cout << "G4TheoFSGenerator: null ptr from transport propagate " << G4endl;
197  throw G4HadronicException(__FILE__, __LINE__, "Null ptr from transport propagate");
198  }
199 } // Uzhi Nov. 2012
200 
201  // Fill particle change
202  unsigned int i;
203  for(i=0; i<theTransportResult->size(); i++)
204  {
205  G4DynamicParticle * aNew =
206  new G4DynamicParticle(theTransportResult->operator[](i)->GetDefinition(),
207  theTransportResult->operator[](i)->GetTotalEnergy(),
208  theTransportResult->operator[](i)->GetMomentum());
209  // @@@ - overkill? G4double newTime = theParticleChange->GetGlobalTime(theTransportResult->operator[](i)->GetFormationTime());
210  theParticleChange->AddSecondary(aNew);
211  delete theTransportResult->operator[](i);
212  }
213 
214  // some garbage collection
215  delete theTransportResult;
216 
217  // Done
218  return theParticleChange;
219 }
G4ReactionProductVector * Propagate(G4KineticTrackVector *theSecondaries, G4V3DNucleus *)
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
virtual G4ReactionProductVector * Propagate(G4KineticTrackVector *theSecondaries, G4V3DNucleus *theNucleus)=0
virtual G4int GetMassNumber()=0
int G4int
Definition: G4Types.hh:78
G4KineticTrackVector * Scatter(G4Nucleus &theNucleus, const G4DynamicParticle &thePrimary)
virtual G4KineticTrackVector * Scatter(const G4Nucleus &theNucleus, const G4DynamicParticle &thePrimary)=0
void SetStatusChange(G4HadFinalStateStatus aS)
void SetPrimaryProjectile(const G4HadProjectile &aPrimary)
virtual G4V3DNucleus * GetWoundedNucleus() const =0
std::vector< G4ReactionProduct * > G4ReactionProductVector
G4IonTable * GetIonTable() const
Hep3Vector vect() const
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
const G4ParticleDefinition * GetDefinition() const
G4double GetIonMass(G4int Z, G4int A, G4int L=0, G4int lvl=0) const
Definition: G4IonTable.cc:1232
G4double GetKineticEnergy() const
static G4Proton * Proton()
Definition: G4Proton.cc:93
virtual G4V3DNucleus * GetProjectileNucleus() const
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
G4double GetFraction(G4Nucleus &theNucleus, const G4DynamicParticle &thePrimary)
const G4LorentzVector & Get4Momentum() const
void SetEnergyChange(G4double anEnergy)
G4double GetPDGMass() const
static G4ParticleTable * GetParticleTable()
virtual G4ReactionProductVector * PropagateNuclNucl(G4KineticTrackVector *theSecondaries, G4V3DNucleus *theNucleus, G4V3DNucleus *theProjectileNucleus)
Hep3Vector unit() const
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
virtual const std::vector< G4Nucleon > & GetNucleons()=0
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void SetMomentumChange(const G4ThreeVector &aV)
void AddSecondary(G4DynamicParticle *aP)
std::pair< G4double, G4double > G4TheoFSGenerator::GetEnergyMomentumCheckLevels ( ) const
virtual

Reimplemented from G4HadronicInteraction.

Definition at line 221 of file G4TheoFSGenerator.cc.

References DBL_MAX, and G4VHighEnergyGenerator::GetEnergyMomentumCheckLevels().

222 {
223  if ( theHighEnergyGenerator ) {
224  return theHighEnergyGenerator->GetEnergyMomentumCheckLevels();
225  } else {
226  return std::pair<G4double, G4double>(DBL_MAX, DBL_MAX);
227  }
228 }
std::pair< G4double, G4double > GetEnergyMomentumCheckLevels() const
#define DBL_MAX
Definition: templates.hh:83
void G4TheoFSGenerator::ModelDescription ( std::ostream &  outFile) const
virtual

Reimplemented from G4HadronicInteraction.

Definition at line 53 of file G4TheoFSGenerator.cc.

References G4VHighEnergyGenerator::GetModelName(), G4HadronicInteraction::GetModelName(), and G4VHighEnergyGenerator::ModelDescription().

54 {
55  outFile << GetModelName() <<" consists of a " << theHighEnergyGenerator->GetModelName()
56  << " string model and of "
57  << ".\n"
58  << "The string model simulates the interaction of\n"
59  << "an incident hadron with a nucleus, forming \n"
60  << "excited strings, decays these strings into hadrons,\n"
61  << "and leaves an excited nucleus.\n"
62  << "The string model:\n";
63  theHighEnergyGenerator->ModelDescription(outFile);
64 //theTransport->IntraNuclearTransportDescription(outFile)
65 }
std::ofstream outFile
Definition: GammaRayTel.cc:68
const G4String & GetModelName() const
virtual G4String GetModelName() const
virtual void ModelDescription(std::ostream &) const
void G4TheoFSGenerator::SetHighEnergyGenerator ( G4VHighEnergyGenerator *const  value)
inline

Definition at line 103 of file G4TheoFSGenerator.hh.

Referenced by G4BertiniElectroNuclearBuilder::Build(), G4FTFBuilder::BuildModel(), G4QGSBuilder::BuildModel(), B03PhysicsList::ConstructHad(), DMXPhysicsList::ConstructHad(), exrdmPhysListHadron::ConstructProcess(), GammaRayTelIonPhysics::ConstructProcess(), GammaRayTelHadronPhysics::ConstructProcess(), G4ElectroVDNuclearModel::G4ElectroVDNuclearModel(), G4FTFBinaryKaonBuilder::G4FTFBinaryKaonBuilder(), G4FTFBinaryNeutronBuilder::G4FTFBinaryNeutronBuilder(), G4FTFBinaryPiKBuilder::G4FTFBinaryPiKBuilder(), G4FTFBinaryPionBuilder::G4FTFBinaryPionBuilder(), G4FTFBinaryProtonBuilder::G4FTFBinaryProtonBuilder(), G4FTFPAntiBarionBuilder::G4FTFPAntiBarionBuilder(), G4FTFPKaonBuilder::G4FTFPKaonBuilder(), G4FTFPNeutronBuilder::G4FTFPNeutronBuilder(), G4FTFPPiKBuilder::G4FTFPPiKBuilder(), G4FTFPPionBuilder::G4FTFPPionBuilder(), G4FTFPProtonBuilder::G4FTFPProtonBuilder(), G4HadronicAbsorptionFritiof::G4HadronicAbsorptionFritiof(), G4HyperonFTFPBuilder::G4HyperonFTFPBuilder(), G4MuonVDNuclearModel::G4MuonVDNuclearModel(), G4QGSBinaryKaonBuilder::G4QGSBinaryKaonBuilder(), G4QGSBinaryNeutronBuilder::G4QGSBinaryNeutronBuilder(), G4QGSBinaryPiKBuilder::G4QGSBinaryPiKBuilder(), G4QGSBinaryPionBuilder::G4QGSBinaryPionBuilder(), G4QGSBinaryProtonBuilder::G4QGSBinaryProtonBuilder(), G4QGSPNeutronBuilder::G4QGSPNeutronBuilder(), G4QGSPPiKBuilder::G4QGSPPiKBuilder(), G4QGSPPionBuilder::G4QGSPPionBuilder(), G4QGSPProtonBuilder::G4QGSPProtonBuilder(), and TheoModelFactory< C, S, F >::New().

104 {
105  theHighEnergyGenerator= value;
106 }
const XML_Char int const XML_Char * value
void G4TheoFSGenerator::SetQuasiElasticChannel ( G4QuasiElasticChannel *const  value)
inline
void G4TheoFSGenerator::SetTransport ( G4VIntraNuclearTransportModel *const  value)
inline

Definition at line 93 of file G4TheoFSGenerator.hh.

Referenced by G4BertiniElectroNuclearBuilder::Build(), G4FTFBuilder::BuildModel(), G4QGSBuilder::BuildModel(), B03PhysicsList::ConstructHad(), DMXPhysicsList::ConstructHad(), exrdmPhysListHadron::ConstructProcess(), GammaRayTelIonPhysics::ConstructProcess(), GammaRayTelHadronPhysics::ConstructProcess(), G4ElectroVDNuclearModel::G4ElectroVDNuclearModel(), G4FTFBinaryKaonBuilder::G4FTFBinaryKaonBuilder(), G4FTFBinaryNeutronBuilder::G4FTFBinaryNeutronBuilder(), G4FTFBinaryPiKBuilder::G4FTFBinaryPiKBuilder(), G4FTFBinaryPionBuilder::G4FTFBinaryPionBuilder(), G4FTFBinaryProtonBuilder::G4FTFBinaryProtonBuilder(), G4FTFPAntiBarionBuilder::G4FTFPAntiBarionBuilder(), G4FTFPKaonBuilder::G4FTFPKaonBuilder(), G4FTFPNeutronBuilder::G4FTFPNeutronBuilder(), G4FTFPPiKBuilder::G4FTFPPiKBuilder(), G4FTFPPionBuilder::G4FTFPPionBuilder(), G4FTFPProtonBuilder::G4FTFPProtonBuilder(), G4HadronicAbsorptionFritiof::G4HadronicAbsorptionFritiof(), G4HyperonFTFPBuilder::G4HyperonFTFPBuilder(), G4MuonVDNuclearModel::G4MuonVDNuclearModel(), G4QGSBinaryKaonBuilder::G4QGSBinaryKaonBuilder(), G4QGSBinaryNeutronBuilder::G4QGSBinaryNeutronBuilder(), G4QGSBinaryPiKBuilder::G4QGSBinaryPiKBuilder(), G4QGSBinaryPionBuilder::G4QGSBinaryPionBuilder(), G4QGSBinaryProtonBuilder::G4QGSBinaryProtonBuilder(), G4QGSPNeutronBuilder::G4QGSPNeutronBuilder(), G4QGSPPiKBuilder::G4QGSPPiKBuilder(), G4QGSPPionBuilder::G4QGSPPionBuilder(), G4QGSPProtonBuilder::G4QGSPProtonBuilder(), and TheoModelFactory< C, S, F >::New().

94 {
95  theTransport = value;
96 }
const XML_Char int const XML_Char * value

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