Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4BinaryLightIonReaction Class Reference

#include <G4BinaryLightIonReaction.hh>

+ Inheritance diagram for G4BinaryLightIonReaction:

Public Member Functions

 G4BinaryLightIonReaction (G4VPreCompoundModel *ptr=0)
 
virtual ~G4BinaryLightIonReaction ()
 
G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &theNucleus)
 
void SetPrecompound (G4VPreCompoundModel *ptr)
 
void SetDeExcitation (G4ExcitationHandler *ptr)
 
virtual void ModelDescription (std::ostream &) const
 
- Public Member Functions inherited from G4HadronicInteraction
 G4HadronicInteraction (const G4String &modelName="HadronicModel")
 
virtual ~G4HadronicInteraction ()
 
virtual G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)=0
 
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
 
G4int GetVerboseLevel () const
 
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, G4doubleGetFatalEnergyCheckLevels () const
 
virtual std::pair< G4double, G4doubleGetEnergyMomentumCheckLevels () const
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
virtual void ModelDescription (std::ostream &outFile) const
 

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 34 of file G4BinaryLightIonReaction.hh.

Constructor & Destructor Documentation

◆ G4BinaryLightIonReaction()

G4BinaryLightIonReaction::G4BinaryLightIonReaction ( G4VPreCompoundModel ptr = 0)

Definition at line 49 of file G4BinaryLightIonReaction.cc.

50: G4HadronicInteraction("Binary Light Ion Cascade"),
51 theProjectileFragmentation(ptr),
52 pA(0),pZ(0), tA(0),tZ(0),spectatorA(0),spectatorZ(0),
53 projectile3dNucleus(0),target3dNucleus(0)
54{
55 if(!ptr) {
58 G4VPreCompoundModel* pre = static_cast<G4VPreCompoundModel*>(p);
59 if(!pre) { pre = new G4PreCompoundModel(); }
60 theProjectileFragmentation = pre;
61 }
62 theModel = new G4BinaryCascade(theProjectileFragmentation);
63 theHandler = theProjectileFragmentation->GetExcitationHandler();
64
65 debug_G4BinaryLightIonReactionResults=getenv("debug_G4BinaryLightIonReactionResults")!=0;
66}
G4HadronicInteraction * FindModel(const G4String &name)
static G4HadronicInteractionRegistry * Instance()
G4ExcitationHandler * GetExcitationHandler() const

◆ ~G4BinaryLightIonReaction()

G4BinaryLightIonReaction::~G4BinaryLightIonReaction ( )
virtual

Definition at line 68 of file G4BinaryLightIonReaction.cc.

69{}

Member Function Documentation

◆ ApplyYourself()

G4HadFinalState * G4BinaryLightIonReaction::ApplyYourself ( const G4HadProjectile aTrack,
G4Nucleus theNucleus 
)
virtual

Implements G4HadronicInteraction.

Definition at line 86 of file G4BinaryLightIonReaction.cc.

88{
89 static G4int eventcounter=0;
90 eventcounter++;
91 if(getenv("BLICDEBUG") ) G4cerr << " ######### Binary Light Ion Reaction number starts ######### "<<eventcounter<<G4endl;
92 G4ping debug("debug_G4BinaryLightIonReaction");
93 pA=aTrack.GetDefinition()->GetBaryonNumber();
94 pZ=G4lrint(aTrack.GetDefinition()->GetPDGCharge()/eplus);
95 tA=targetNucleus.GetA_asInt();
96 tZ=targetNucleus.GetZ_asInt();
97
98 G4LorentzVector mom(aTrack.Get4Momentum());
99 //G4cout << "proj mom : " << mom << G4endl;
100 G4LorentzRotation toBreit(mom.boostVector());
101
102 G4bool swapped=SetLighterAsProjectile(mom, toBreit);
103 //G4cout << "after swap, swapped? / mom " << swapped << " / " << mom <<G4endl;
104 G4ReactionProductVector * result = 0;
105 G4ReactionProductVector * cascaders=0; //new G4ReactionProductVector;
106// G4double m_nucl(0); // to check energy balance
107
108
109 // G4double m1=G4ParticleTable::GetParticleTable()->GetIonTable()->GetIonMass(pZ,pA);
110 // G4cout << "Entering the decision point "
111 // << (mom.t()-mom.mag())/pA << " "
112 // << pA<<" "<< pZ<<" "
113 // << tA<<" "<< tZ<<G4endl
114 // << " "<<mom.t()-mom.mag()<<" "
115 // << mom.t()- m1<<G4endl;
116 if( (mom.t()-mom.mag())/pA < 50*MeV )
117 {
118 // G4cout << "Using pre-compound only, E= "<<mom.t()-mom.mag()<<G4endl;
119 // m_nucl = mom.mag();
120 cascaders=FuseNucleiAndPrompound(mom);
121 }
122 else
123 {
124 result=Interact(mom,toBreit);
125
126 if(! result )
127 {
128 {
129 // abort!!
130
131 G4cerr << "G4BinaryLightIonReaction no final state for: " << G4endl;
132 G4cerr << " Primary " << aTrack.GetDefinition()
133 << ", (A,Z)=(" << aTrack.GetDefinition()->GetBaryonNumber()
134 << "," << aTrack.GetDefinition()->GetPDGCharge()/eplus << ") "
135 << ", kinetic energy " << aTrack.GetKineticEnergy()
136 << G4endl;
137 G4cerr << " Target nucleus (A,Z)=("
138 << (swapped?pA:tA) << ","
139 << (swapped?pZ:tZ) << ")" << G4endl;
140 G4cerr << " if frequent, please submit above information as bug report"
141 << G4endl << G4endl;
142
143 theResult.Clear();
144 theResult.SetStatusChange(isAlive);
145 theResult.SetEnergyChange(aTrack.GetKineticEnergy());
146 theResult.SetMomentumChange(aTrack.Get4Momentum().vect().unit());
147 return &theResult;
148
149 }
150 }
151
152 // Calculate excitation energy,
153 G4double theStatisticalExEnergy = GetProjectileExcitation();
154
155
156 pInitialState = mom;
157 //G4cout << "pInitialState from aTrack : " << pInitialState;
158 pInitialState.setT(pInitialState.getT() +
160 //G4cout << "target nucleus added : " << pInitialState << G4endl;
161
162 delete target3dNucleus;target3dNucleus=0;
163 delete projectile3dNucleus;projectile3dNucleus=0;
164
166
167 cascaders = new G4ReactionProductVector;
168
169 G4LorentzVector pspectators=SortResult(result,spectators,cascaders);
170
171 // pFinalState=std::accumulate(cascaders->begin(),cascaders->end(),pFinalState,ReactionProduct4Mom);
172 std::vector<G4ReactionProduct *>::iterator iter;
173
174 //G4cout << "pInitialState, pFinalState / pspectators"<< pInitialState << " / " << pFinalState << " / " << pspectators << G4endl;
175 // if ( spectA-spectatorA !=0 || spectZ-spectatorZ !=0)
176 // {
177 // G4cout << "spect Nucl != spectators: nucl a,z; spect a,z" <<
178 // spectatorA <<" "<< spectatorZ <<" ; " << spectA <<" "<< spectZ << G4endl;
179 // }
180 delete result;
181 result=0;
182 G4LorentzVector momentum(pInitialState-pFinalState);
183 G4int loopcount(0);
184 //G4cout << "momentum, pspectators : " << momentum << " / " << pspectators << G4endl;
185 while (std::abs(momentum.e()-pspectators.e()) > 10*MeV)
186 {
187 G4LorentzVector pCorrect(pInitialState-pspectators);
188 //G4cout << "BIC nonconservation? (pInitialState-pFinalState) / spectators :" << momentum << " / " << pspectators << "pCorrect "<< pCorrect<< G4endl;
189 // Correct outgoing casacde particles.... to have momentum of (initial state - spectators)
190 G4bool EnergyIsCorrect=EnergyAndMomentumCorrector(cascaders, pCorrect);
191 if ( ! EnergyIsCorrect && debug_G4BinaryLightIonReactionResults)
192 {
193 G4cout << "Warning - G4BinaryLightIonReaction E/P correction for cascaders failed" << G4endl;
194 }
195 pFinalState=G4LorentzVector(0,0,0,0);
196 unsigned int i;
197 for(i=0; i<cascaders->size(); i++)
198 {
199 pFinalState += G4LorentzVector( (*cascaders)[i]->GetMomentum(), (*cascaders)[i]->GetTotalEnergy() );
200 }
201 momentum=pInitialState-pFinalState;
202 if (++loopcount > 10 )
203 {
204 if ( momentum.vect().mag() - momentum.e()> 10*keV )
205 {
206 G4cerr << "G4BinaryLightIonReaction.cc: Cannot correct 4-momentum of cascade particles" << G4endl;
207 throw G4HadronicException(__FILE__, __LINE__, "G4BinaryCasacde::ApplyCollision()");
208 } else {
209 break;
210 }
211
212 }
213 }
214
215 if (spectatorA > 0 )
216 {
217 // check spectator momentum
218 if ( momentum.vect().mag() - momentum.e()> 10*keV )
219 {
220
221 G4ReactionProductVector::iterator ispectator;
222 for (ispectator=spectators->begin();ispectator!=spectators->end();ispectator++)
223 {
224 delete *ispectator;
225 }
226 delete spectators;
227
228 G4cout << "G4BinaryLightIonReaction.cc: mom check: " << momentum
229 << " 3.mag "<< momentum.vect().mag() << G4endl
230 << " .. pInitialState/pFinalState/spectators " << pInitialState <<" "
231 << pFinalState << " " << pspectators << G4endl
232 << " .. A,Z " << spectatorA <<" "<< spectatorZ << G4endl;
233 G4cout << "G4BinaryLightIonReaction invalid final state for: " << G4endl;
234 G4cout << " Primary " << aTrack.GetDefinition()
235 << ", (A,Z)=(" << aTrack.GetDefinition()->GetBaryonNumber()
236 << "," << aTrack.GetDefinition()->GetPDGCharge()/eplus << ") "
237 << ", kinetic energy " << aTrack.GetKineticEnergy()
238 << G4endl;
239 G4cout << " Target nucleus (A,Z)=(" << targetNucleus.GetA_asInt()
240 << "," << targetNucleus.GetZ_asInt() << ")" << G4endl;
241 G4cout << " if frequent, please submit above information as bug report"
242 << G4endl << G4endl;
243
244 theResult.Clear();
245 theResult.SetStatusChange(isAlive);
246 theResult.SetEnergyChange(aTrack.GetKineticEnergy());
247 theResult.SetMomentumChange(aTrack.Get4Momentum().vect().unit());
248 return &theResult;
249 }
250
251
252 DeExciteSpectatorNucleus(spectators, cascaders, theStatisticalExEnergy, momentum);
253 }
254 }
255 // Rotate to lab
257 toZ.rotateZ(-1*mom.phi());
258 toZ.rotateY(-1*mom.theta());
259 G4LorentzRotation toLab(toZ.inverse());
260
261 // Fill the particle change, while rotating. Boost from projectile breit-frame in case we swapped.
262 // theResult.Clear();
263 theResult.Clear();
264 theResult.SetStatusChange(stopAndKill);
265 G4double Etot(0);
266 size_t i=0;
267 for(i=0; i<cascaders->size(); i++)
268 {
269 if((*cascaders)[i]->GetNewlyAdded())
270 {
271 G4DynamicParticle * aNew =
272 new G4DynamicParticle((*cascaders)[i]->GetDefinition(),
273 (*cascaders)[i]->GetTotalEnergy(),
274 (*cascaders)[i]->GetMomentum() );
275 G4LorentzVector tmp = aNew->Get4Momentum();
276 if(swapped)
277 {
278 tmp*=toBreit.inverse();
279 tmp.setVect(-tmp.vect());
280 }
281 tmp *= toLab;
282 aNew->Set4Momentum(tmp);
283 //G4cout << "result[" << i << "], 4vect: " << tmp << G4endl;
284 theResult.AddSecondary(aNew);
285 Etot += tmp.e();
286 // G4cout << "LIBIC: Secondary " << aNew->GetDefinition()->GetParticleName()
287 // <<" "<< aNew->GetMomentum()
288 // <<" "<< aNew->GetTotalEnergy()
289 // << G4endl;
290 }
291 delete (*cascaders)[i];
292 }
293 delete cascaders;
294
295#ifdef debug_BLIR_result
296 G4cout << "Result analysis, secondaries" << theResult.GetNumberOfSecondaries() << G4endl;
297 G4cout << " Energy conservation initial/primary/nucleus/final/delta(init-final) "
298 << aTrack.GetTotalEnergy() + m_nucl << aTrack.GetTotalEnergy() << m_nucl <<Etot
299 << aTrack.GetTotalEnergy() + m_nucl - Etot;
300#endif
301
302 if(getenv("BLICDEBUG") ) G4cerr << " ######### Binary Light Ion Reaction number ends ######### "<<eventcounter<<G4endl;
303
304 return &theResult;
305}
@ isAlive
@ stopAndKill
CLHEP::HepLorentzVector G4LorentzVector
std::vector< G4ReactionProduct * > G4ReactionProductVector
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cerr
G4DLLIMPORT std::ostream G4cout
Hep3Vector unit() const
HepLorentzRotation & rotateY(double delta)
HepLorentzRotation & rotateZ(double delta)
HepLorentzRotation inverse() const
Hep3Vector vect() const
void setVect(const Hep3Vector &)
G4LorentzVector Get4Momentum() const
void Set4Momentum(const G4LorentzVector &momentum)
void SetStatusChange(G4HadFinalStateStatus aS)
G4int GetNumberOfSecondaries() const
void AddSecondary(G4DynamicParticle *aP)
void SetEnergyChange(G4double anEnergy)
void SetMomentumChange(const G4ThreeVector &aV)
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
G4double GetTotalEnergy() const
G4double GetIonMass(G4int Z, G4int A, G4int L=0) const
!! Only ground states are supported now
Definition: G4IonTable.cc:774
G4double GetPDGCharge() const
static G4ParticleTable * GetParticleTable()
G4IonTable * GetIonTable()
Definition: G4ping.hh:34
int G4lrint(double ad)
Definition: templates.hh:163

◆ ModelDescription()

void G4BinaryLightIonReaction::ModelDescription ( std::ostream &  outFile) const
virtual

Reimplemented from G4HadronicInteraction.

Definition at line 71 of file G4BinaryLightIonReaction.cc.

72{
73 outFile << "G4Binary Light Ion Cascade is an intra-nuclear cascade model\n"
74 << "using G4BinaryCasacde to model the interaction of a light\n"
75 << "nucleus with a nucleus.\n"
76 << "The lighter of the two nuclei is treated like a set of projectiles\n"
77 << "which are transported simultanously through the heavier nucleus.\n";
78}

◆ SetDeExcitation()

void G4BinaryLightIonReaction::SetDeExcitation ( G4ExcitationHandler ptr)
inline

Definition at line 74 of file G4BinaryLightIonReaction.hh.

75{
76 theProjectileFragmentation->SetExcitationHandler(ptr);
77 theHandler = ptr;
78}
void SetExcitationHandler(G4ExcitationHandler *ptr)

◆ SetPrecompound()

void G4BinaryLightIonReaction::SetPrecompound ( G4VPreCompoundModel ptr)
inline

Definition at line 69 of file G4BinaryLightIonReaction.hh.

70{
71 if(ptr) { theProjectileFragmentation = ptr; }
72 theHandler = theProjectileFragmentation->GetExcitationHandler();
73}

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