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

#include <G4RPGFragmentation.hh>

+ Inheritance diagram for G4RPGFragmentation:

Public Member Functions

 G4RPGFragmentation ()
 
void FragmentationIntegral (G4double, G4double, G4double, G4double)
 
G4bool ReactionStage (const G4HadProjectile *, G4ReactionProduct &, G4bool &, const G4DynamicParticle *, G4ReactionProduct &, G4bool &, const G4Nucleus &, G4ReactionProduct &, G4FastVector< G4ReactionProduct, 256 > &, G4int &, G4bool, G4ReactionProduct &)
 
- Public Member Functions inherited from G4RPGReaction
 G4RPGReaction ()
 
virtual ~G4RPGReaction ()
 
G4bool ReactionStage (const G4HadProjectile *, G4ReactionProduct &, G4bool &, const G4DynamicParticle *, G4ReactionProduct &, G4bool &, const G4Nucleus &, G4ReactionProduct &, G4FastVector< G4ReactionProduct, 256 > &, G4int &, G4bool, G4ReactionProduct &)
 
void AddBlackTrackParticles (const G4double, const G4int, const G4double, const G4int, const G4ReactionProduct &, G4int, G4int, const G4Nucleus &, G4FastVector< G4ReactionProduct, 256 > &, G4int &)
 
G4double GenerateNBodyEvent (const G4double totalEnergy, const G4bool constantCrossSection, G4FastVector< G4ReactionProduct, 256 > &vec, G4int &vecLen)
 
G4double GenerateNBodyEventT (const G4double totalEnergy, const G4bool constantCrossSection, std::vector< G4ReactionProduct * > &list)
 
void NuclearReaction (G4FastVector< G4ReactionProduct, 4 > &vec, G4int &vecLen, const G4HadProjectile *originalIncident, const G4Nucleus &aNucleus, const G4double theAtomicMass, const G4double *massVec)
 

Additional Inherited Members

- Protected Member Functions inherited from G4RPGReaction
void Rotate (const G4double numberofFinalStateNucleons, const G4ThreeVector &temp, const G4ReactionProduct &modifiedOriginal, const G4HadProjectile *originalIncident, const G4Nucleus &targetNucleus, G4ReactionProduct &currentParticle, G4ReactionProduct &targetParticle, G4FastVector< G4ReactionProduct, 256 > &vec, G4int &vecLen)
 
void Defs1 (const G4ReactionProduct &modifiedOriginal, G4ReactionProduct &currentParticle, G4ReactionProduct &targetParticle, G4FastVector< G4ReactionProduct, 256 > &vec, G4int &vecLen)
 
std::pair< G4int, G4intGetFinalStateNucleons (const G4DynamicParticle *originalTarget, const G4FastVector< G4ReactionProduct, 256 > &vec, const G4int &vecLen)
 
void MomentumCheck (const G4ReactionProduct &modifiedOriginal, G4ReactionProduct &currentParticle, G4ReactionProduct &targetParticle, G4FastVector< G4ReactionProduct, 256 > &vec, G4int &vecLen)
 
G4double normal ()
 
G4ThreeVector Isotropic (const G4double &)
 

Detailed Description

Definition at line 48 of file G4RPGFragmentation.hh.

Constructor & Destructor Documentation

◆ G4RPGFragmentation()

G4RPGFragmentation::G4RPGFragmentation ( )

Definition at line 42 of file G4RPGFragmentation.cc.

44{
45 for (G4int i = 0; i < 20; i++) dndl[i] = 0.0;
46}
int G4int
Definition: G4Types.hh:85

Member Function Documentation

◆ FragmentationIntegral()

void G4RPGFragmentation::FragmentationIntegral ( G4double  pt,
G4double  et,
G4double  parMass,
G4double  secMass 
)

Definition at line 49 of file G4RPGFragmentation.cc.

51{
52 pt = std::max( 0.001, pt );
53 G4double dx = 1./(19.*pt);
54 G4double x;
55 G4double term1;
56 G4double term2;
57
58 for (G4int i = 1; i < 20; i++) {
59 x = (G4double(i) - 0.5)*dx;
60 term1 = 1. + parMass*parMass*x*x;
61 term2 = pt*x*et*pt*x*et + pt*pt + secMass*secMass;
62 dndl[i] = dx / std::sqrt( term1*term1*term1*term2 )
63 + dndl[i-1];
64 }
65}
double G4double
Definition: G4Types.hh:83

Referenced by ReactionStage().

◆ ReactionStage()

G4bool G4RPGFragmentation::ReactionStage ( const G4HadProjectile originalIncident,
G4ReactionProduct modifiedOriginal,
G4bool incidentHasChanged,
const G4DynamicParticle originalTarget,
G4ReactionProduct targetParticle,
G4bool targetHasChanged,
const G4Nucleus targetNucleus,
G4ReactionProduct currentParticle,
G4FastVector< G4ReactionProduct, 256 > &  vec,
G4int vecLen,
G4bool  leadFlag,
G4ReactionProduct leadingStrangeParticle 
)

Definition at line 68 of file G4RPGFragmentation.cc.

81{
82 //
83 // Based on H. Fesefeldt's original FORTRAN code GENXPT
84 //
85 // Generation of x- and pT- values for incident, target, and all secondary
86 // particles using a simple single variable description E D3S/DP3= F(Q)
87 // with Q^2 = (M*X)^2 + PT^2. Final state kinematics are produced by an
88 // FF-type iterative cascade method
89 //
90 // Internal units are GeV
91 //
92
93 // Protection in case no secondary has been created. In that case use
94 // two-body scattering
95 //
96 if (vecLen == 0) return false;
97
103
104 G4int i, l;
105 G4bool veryForward = false;
106
107 const G4double ekOriginal = modifiedOriginal.GetKineticEnergy()/GeV;
108 const G4double etOriginal = modifiedOriginal.GetTotalEnergy()/GeV;
109 const G4double mOriginal = modifiedOriginal.GetMass()/GeV;
110 const G4double pOriginal = modifiedOriginal.GetMomentum().mag()/GeV;
111 G4double targetMass = targetParticle.GetDefinition()->GetPDGMass()/GeV;
112 G4double centerofmassEnergy = std::sqrt( mOriginal*mOriginal +
113 targetMass*targetMass +
114 2.0*targetMass*etOriginal ); // GeV
115 G4double currentMass = currentParticle.GetMass()/GeV;
116 targetMass = targetParticle.GetMass()/GeV;
117
118 // Randomize the order of the secondary particles.
119 // Note that the current and target particles are not affected.
120
121 for (i=0; i<vecLen; ++i) {
122 G4int itemp = G4int( G4UniformRand()*vecLen );
123 G4ReactionProduct pTemp = *vec[itemp];
124 *vec[itemp] = *vec[i];
125 *vec[i] = pTemp;
126 }
127
128 if (currentMass == 0.0 && targetMass == 0.0) {
129 // Target and projectile have annihilated. Replace them with the first
130 // two secondaries in the list. Current particle KE is maintained.
131
132 G4double ek = currentParticle.GetKineticEnergy();
133 G4ThreeVector mom = currentParticle.GetMomentum();
134 currentParticle = *vec[0];
135 currentParticle.SetSide(1);
136 targetParticle = *vec[1];
137 targetParticle.SetSide(-1);
138 for( i=0; i<(vecLen-2); ++i )*vec[i] = *vec[i+2];
139 G4ReactionProduct *temp = vec[vecLen-1];
140 delete temp;
141 temp = vec[vecLen-2];
142 delete temp;
143 vecLen -= 2;
144 currentMass = currentParticle.GetMass()/GeV;
145 targetMass = targetParticle.GetMass()/GeV;
146 incidentHasChanged = true;
147 targetHasChanged = true;
148 currentParticle.SetKineticEnergy( ek );
149 currentParticle.SetMomentum(mom);
150 veryForward = true;
151 }
152 const G4double atomicWeight = targetNucleus.GetA_asInt();
153 const G4double atomicNumber = targetNucleus.GetZ_asInt();
154 const G4double protonMass = aProton->GetPDGMass();
155
156 if (originalIncident->GetDefinition()->GetParticleSubType() == "kaon"
157 && G4UniformRand() >= 0.7) {
158 G4ReactionProduct temp = currentParticle;
159 currentParticle = targetParticle;
160 targetParticle = temp;
161 incidentHasChanged = true;
162 targetHasChanged = true;
163 currentMass = currentParticle.GetMass()/GeV;
164 targetMass = targetParticle.GetMass()/GeV;
165 }
166 const G4double afc = std::min( 0.75,
167 0.312+0.200*G4Log(G4Log(centerofmassEnergy*centerofmassEnergy))+
168 std::pow(centerofmassEnergy*centerofmassEnergy,1.5)/6000.0 );
169
170 G4double freeEnergy = centerofmassEnergy-currentMass-targetMass;
171 G4double forwardEnergy = freeEnergy/2.;
172 G4int forwardCount = 1; // number of particles in forward hemisphere
173
174 G4double backwardEnergy = freeEnergy/2.;
175 G4int backwardCount = 1; // number of particles in backward hemisphere
176
177 if(veryForward) {
178 if(currentParticle.GetSide()==-1)
179 {
180 forwardEnergy += currentMass;
181 forwardCount --;
182 backwardEnergy -= currentMass;
183 backwardCount ++;
184 }
185 if(targetParticle.GetSide()!=-1)
186 {
187 backwardEnergy += targetMass;
188 backwardCount --;
189 forwardEnergy -= targetMass;
190 forwardCount ++;
191 }
192 }
193
194 for (i=0; i<vecLen; ++i) {
195 if( vec[i]->GetSide() == -1 )
196 {
197 ++backwardCount;
198 backwardEnergy -= vec[i]->GetMass()/GeV;
199 } else {
200 ++forwardCount;
201 forwardEnergy -= vec[i]->GetMass()/GeV;
202 }
203 }
204
205 // Check that sum of forward particle masses does not exceed forwardEnergy,
206 // and similarly for backward particles. If so, move particles from one
207 // hemisphere to another.
208
209 if (backwardEnergy < 0.0) {
210 for (i = 0; i < vecLen; ++i) {
211 if (vec[i]->GetSide() == -1) {
212 backwardEnergy += vec[i]->GetMass()/GeV;
213 --backwardCount;
214 vec[i]->SetSide(1);
215 forwardEnergy -= vec[i]->GetMass()/GeV;
216 ++forwardCount;
217 if (backwardEnergy > 0.0) break;
218 }
219 }
220 }
221
222 if (forwardEnergy < 0.0) {
223 for (i = 0; i < vecLen; ++i) {
224 if (vec[i]->GetSide() == 1) {
225 forwardEnergy += vec[i]->GetMass()/GeV;
226 --forwardCount;
227 vec[i]->SetSide(-1);
228 backwardEnergy -= vec[i]->GetMass()/GeV;
229 ++backwardCount;
230 if (forwardEnergy > 0.0) break;
231 }
232 }
233 }
234
235 // Special cases for reactions near threshold
236
237 // 1. There is only one secondary
238 if (forwardEnergy > 0.0 && backwardEnergy < 0.0) {
239 forwardEnergy += backwardEnergy;
240 backwardEnergy = 0;
241 }
242
243 // 2. Nuclear binding energy is large
244 if (forwardEnergy + backwardEnergy < 0.0) return false;
245
246
247 // forwardEnergy now represents the total energy in the forward reaction
248 // hemisphere which is available for kinetic energy and particle creation.
249 // Similarly for backwardEnergy.
250
251 // Add particles from the intranuclear cascade.
252 // nuclearExcitationCount = number of new secondaries produced by nuclear
253 // excitation
254 // extraCount = number of nucleons within these new secondaries
255 //
256 // Note: eventually have to make sure that enough nucleons are available
257 // in the case of small target nuclei
258
259 G4double xtarg;
260 G4double a13 = G4Pow::GetInstance()->A13(atomicWeight); // A**(1/3)
261 if (centerofmassEnergy < (2.0+G4UniformRand()) )
262 xtarg = afc * (a13-1.0) * (2.0*backwardCount+vecLen+2)/2.0;
263 else
264 xtarg = afc * (a13-1.0) * (2.0*backwardCount);
265
266 if( xtarg <= 0.0 )xtarg = 0.01;
267 G4int nuclearExcitationCount = G4Poisson( xtarg );
268 // To do: try reducing nuclearExcitationCount with increasing energy
269 // to simulate cut-off of cascade
270 if(atomicWeight<1.0001) nuclearExcitationCount = 0;
271 G4int extraNucleonCount = 0;
272 G4double extraNucleonMass = 0.0;
273
274 if (nuclearExcitationCount > 0) {
275 const G4double nucsup[] = { 1.00, 0.7, 0.5, 0.4, 0.35, 0.3 };
276 const G4double psup[] = { 3., 6., 20., 50., 100., 1000. };
277 G4int momentumBin = 0;
278
279 G4int loop = 0;
281 ed << " While count exceeded " << G4endl;
282 while( (momentumBin < 6) &&
283 (modifiedOriginal.GetTotalMomentum()/GeV > psup[momentumBin]) ) { /* Loop checking, 01.09.2015, D.Wright */
284 ++momentumBin;
285 loop++;
286 if (loop > 1000) {
287 G4Exception("G4RPGFragmentation::ReactionStage()", "HAD_RPG_100", JustWarning, ed);
288 break;
289 }
290 }
291
292 momentumBin = std::min( 5, momentumBin );
293
294 // NOTE: in GENXPT, these new particles were given negative codes
295 // here I use NewlyAdded = true instead
296 //
297
298 for (i = 0; i < nuclearExcitationCount; ++i) {
300 if (G4UniformRand() < nucsup[momentumBin]) {
301
302 if (G4UniformRand() > 1.0-atomicNumber/atomicWeight)
303 pVec->SetDefinition( aProton );
304 else
305 pVec->SetDefinition( aNeutron );
306
307 // nucleon comes from nucleus -
308 // do not subtract its mass from backward energy
309 pVec->SetSide( -2 ); // -2 means backside nucleon
310 ++extraNucleonCount;
311 extraNucleonMass += pVec->GetMass()/GeV;
312 // To do: Remove chosen nucleon from target nucleus
313 pVec->SetNewlyAdded( true );
314 vec.SetElement( vecLen++, pVec );
315 ++backwardCount;
316
317 } else {
318
319 G4double ran = G4UniformRand();
320 if( ran < 0.3181 )
321 pVec->SetDefinition( aPiPlus );
322 else if( ran < 0.6819 )
323 pVec->SetDefinition( aPiZero );
324 else
325 pVec->SetDefinition( aPiMinus );
326
327 if (backwardEnergy > pVec->GetMass()/GeV) {
328 backwardEnergy -= pVec->GetMass()/GeV; // pion mass comes from free energy
329 ++backwardCount;
330 pVec->SetSide( -1 ); // backside particle, but not a nucleon
331 pVec->SetNewlyAdded( true );
332 vec.SetElement( vecLen++, pVec );
333 }
334
335 // To do: Change proton to neutron (or vice versa) in target nucleus depending
336 // on pion charge
337 }
338 }
339 }
340
341 // Define initial state vectors for Lorentz transformations
342 // The pseudoParticles have non-standard masses, hence the "pseudo"
343
344 G4ReactionProduct pseudoParticle[8];
345 for (i = 0; i < 8; ++i) pseudoParticle[i].SetZero();
346
347 pseudoParticle[0].SetMass( mOriginal*GeV );
348 pseudoParticle[0].SetMomentum( 0.0, 0.0, pOriginal*GeV );
349 pseudoParticle[0].SetTotalEnergy(
350 std::sqrt( pOriginal*pOriginal + mOriginal*mOriginal )*GeV );
351
352 pseudoParticle[1].SetMass(protonMass); // this could be targetMass
353 pseudoParticle[1].SetTotalEnergy(protonMass);
354
355 pseudoParticle[3].SetMass(protonMass*(1+extraNucleonCount) );
356 pseudoParticle[3].SetTotalEnergy(protonMass*(1+extraNucleonCount) );
357
358 pseudoParticle[2] = pseudoParticle[0] + pseudoParticle[1];
359 pseudoParticle[3] = pseudoParticle[3] + pseudoParticle[0];
360
361 pseudoParticle[0].Lorentz( pseudoParticle[0], pseudoParticle[2] );
362 pseudoParticle[1].Lorentz( pseudoParticle[1], pseudoParticle[2] );
363
364 // Main loop for 4-momentum generation
365 // See Pitha-report (Aachen) for a detailed description of the method
366
367 G4double aspar, pt, et, x, pp, pp1, wgt;
368 G4int innerCounter, outerCounter;
369 G4bool eliminateThisParticle, resetEnergies, constantCrossSection;
370
371 G4double forwardKinetic = 0.0;
372 G4double backwardKinetic = 0.0;
373
374 // Process the secondary particles in reverse order
375 // The incident particle is done after the secondaries
376 // Nucleons, including the target, in the backward hemisphere are also
377 // done later
378
379 G4int backwardNucleonCount = 0; // number of nucleons in backward hemisphere
380 G4double totalEnergy, kineticEnergy, vecMass;
381 G4double phi;
382
383 for (i = vecLen-1; i >= 0; --i) {
384
385 if (vec[i]->GetNewlyAdded()) { // added from intranuclear cascade
386 if (vec[i]->GetSide() == -2) { // its a nucleon
387 if (backwardNucleonCount < 18) {
388 if (vec[i]->GetDefinition()->GetParticleSubType() == "pi") {
389 for (G4int j = 0; j < vecLen; j++) delete vec[j];
390 vecLen = 0;
391 throw G4HadReentrentException(__FILE__, __LINE__,
392 "G4RPGFragmentation::ReactionStage : a pion has been counted as a backward nucleon");
393 }
394 vec[i]->SetSide(-3);
395 ++backwardNucleonCount;
396 continue; // Don't generate momenta for the first 17 backward
397 // cascade nucleons. This gets done by the cluster
398 // model later on.
399 }
400 }
401 }
402
403 // Set pt and phi values, they are changed somewhat in the iteration loop
404 // Set mass parameter for lambda fragmentation model
405
406 vecMass = vec[i]->GetMass()/GeV;
407 G4double ran = -G4Log(1.0-G4UniformRand())/3.5;
408
409 if (vec[i]->GetSide() == -2) { // backward nucleon
410 aspar = 0.20;
411 pt = std::sqrt( std::pow( ran, 1.2 ) );
412
413 } else { // not a backward nucleon
414 if (vec[i]->GetDefinition()->GetParticleSubType() == "pi") {
415 aspar = 0.75;
416 pt = std::sqrt( std::pow( ran, 1.7 ) );
417 } else if (vec[i]->GetDefinition()->GetParticleSubType() == "kaon") {
418 aspar = 0.70;
419 pt = std::sqrt( std::pow( ran, 1.7 ) );
420 } else { // vec[i] must be a baryon or ion
421 aspar = 0.65;
422 pt = std::sqrt( std::pow( ran, 1.5 ) );
423 }
424 }
425
426 pt = std::max( 0.001, pt );
427 phi = G4UniformRand()*twopi;
428 vec[i]->SetMomentum( pt*std::cos(phi)*GeV, pt*std::sin(phi)*GeV );
429 if (vec[i]->GetSide() > 0)
430 et = pseudoParticle[0].GetTotalEnergy()/GeV;
431 else
432 et = pseudoParticle[1].GetTotalEnergy()/GeV;
433
434 //
435 // Start of outer iteration loop
436 //
437 outerCounter = 0;
438 eliminateThisParticle = true;
439 resetEnergies = true;
440 dndl[0] = 0.0;
441
442 while (++outerCounter < 3) { /* Loop checking, 01.09.2015, D.Wright */
443
444 FragmentationIntegral(pt, et, aspar, vecMass);
445 innerCounter = 0;
446 vec[i]->SetMomentum( pt*std::cos(phi)*GeV, pt*std::sin(phi)*GeV );
447
448 // Start of inner iteration loop
449
450 while (++innerCounter < 7) { /* Loop checking, 01.09.2015, D.Wright */
451
452 ran = G4UniformRand()*dndl[19];
453 l = 1;
454
455 G4int loop = 0;
457 ed << " While count exceeded " << G4endl;
458 while( ( ran > dndl[l] ) && ( l < 19 ) ) { /* Loop checking, 01.09.2015, D.Wright */
459 l++;
460 loop++;
461 if (loop > 1000) {
462 G4Exception("G4RPGFragmentation::ReactionStage()", "HAD_RPG_100", JustWarning, ed);
463 break;
464 }
465 }
466
467 x = (G4double(l-1) + G4UniformRand())/19.;
468 if (vec[i]->GetSide() < 0) x *= -1.;
469 vec[i]->SetMomentum( x*et*GeV ); // set the z-momentum
470 totalEnergy = std::sqrt( x*et*x*et + pt*pt + vecMass*vecMass );
471 vec[i]->SetTotalEnergy( totalEnergy*GeV );
472 kineticEnergy = vec[i]->GetKineticEnergy()/GeV;
473
474 if (vec[i]->GetSide() > 0) { // forward side
475 if( (forwardKinetic+kineticEnergy) < 0.95*forwardEnergy ) {
476 // Leave at least 5% of the forward free energy for the projectile
477
478 pseudoParticle[4] = pseudoParticle[4] + (*vec[i]);
479 forwardKinetic += kineticEnergy;
480 outerCounter = 2; // leave outer loop
481 eliminateThisParticle = false; // don't eliminate this particle
482 resetEnergies = false;
483 break; // leave inner loop
484 }
485 if( innerCounter > 5 )break; // leave inner loop
486 if( backwardEnergy >= vecMass ) // switch sides
487 {
488 vec[i]->SetSide(-1);
489 forwardEnergy += vecMass;
490 backwardEnergy -= vecMass;
491 ++backwardCount;
492 }
493 } else { // backward side
494 // if (extraNucleonCount > 19) x = 0.999;
495 // G4double xxx = 0.95+0.05*extraNucleonCount/20.0;
496 // DHW: I think above lines were meant to be as follows:
497 G4double xxx = 0.999;
498 if (extraNucleonCount < 20) xxx = 0.95+0.05*extraNucleonCount/20.0;
499
500 if ((backwardKinetic+kineticEnergy) < xxx*backwardEnergy) {
501 pseudoParticle[5] = pseudoParticle[5] + (*vec[i]);
502 backwardKinetic += kineticEnergy;
503 outerCounter = 2; // leave outer loop
504 eliminateThisParticle = false; // don't eliminate this particle
505 resetEnergies = false;
506 break; // leave inner loop
507 }
508 if (innerCounter > 5) break; // leave inner loop
509 if (forwardEnergy >= vecMass) { // switch sides
510 vec[i]->SetSide(1);
511 forwardEnergy -= vecMass;
512 backwardEnergy += vecMass;
513 backwardCount--;
514 }
515 }
516 G4ThreeVector momentum = vec[i]->GetMomentum();
517 vec[i]->SetMomentum( momentum.x() * 0.9, momentum.y() * 0.9 );
518 pt *= 0.9;
519 dndl[19] *= 0.9;
520 } // closes inner loop
521
522 // If we get here, the inner loop has been done 6 times.
523 // If necessary, reduce energies of the previously done particles if
524 // they are lighter than protons or are in the forward hemisphere.
525 // Then continue with outer loop.
526
527 if (resetEnergies)
528 ReduceEnergiesOfSecondaries(i+1, forwardKinetic, backwardKinetic,
529 vec, vecLen,
530 pseudoParticle[4], pseudoParticle[5],
531 pt);
532
533 } // closes outer loop
534
535 if (eliminateThisParticle && vec[i]->GetMayBeKilled()) {
536 // not enough energy, eliminate this particle
537
538 if (vec[i]->GetSide() > 0) {
539 --forwardCount;
540 forwardEnergy += vecMass;
541 } else {
542 --backwardCount;
543 if (vec[i]->GetSide() == -2) {
544 --extraNucleonCount;
545 extraNucleonMass -= vecMass;
546 } else {
547 backwardEnergy += vecMass;
548 }
549 }
550
551 for( G4int j=i; j<(vecLen-1); ++j )*vec[j] = *vec[j+1]; // shift up
552 G4ReactionProduct* temp = vec[vecLen-1];
553 delete temp;
554 // To do: modify target nucleus according to particle eliminated
555
556 if( --vecLen == 0 ){
557 G4cout << " FALSE RETURN DUE TO ENERGY BALANCE " << G4endl;
558 return false;
559 } // all the secondaries have been eliminated
560 }
561 } // closes main loop over secondaries
562
563 // Re-balance forward and backward energy if possible and if necessary
564
565 G4double forwardKEDiff = forwardEnergy - forwardKinetic;
566 G4double backwardKEDiff = backwardEnergy - backwardKinetic;
567
568 if (forwardKEDiff < 0.0 || backwardKEDiff < 0.0) {
569 ReduceEnergiesOfSecondaries(0, forwardKinetic, backwardKinetic,
570 vec, vecLen,
571 pseudoParticle[4], pseudoParticle[5],
572 pt);
573
574 forwardKEDiff = forwardEnergy - forwardKinetic;
575 backwardKEDiff = backwardEnergy - backwardKinetic;
576 if (backwardKEDiff < 0.0) {
577 if (forwardKEDiff + backwardKEDiff > 0.0) {
578 backwardEnergy = backwardKinetic;
579 forwardEnergy += backwardKEDiff;
580 forwardKEDiff = forwardEnergy - forwardKinetic;
581 backwardKEDiff = 0.0;
582 } else {
583 G4cout << " False return due to insufficient backward energy " << G4endl;
584 return false;
585 }
586 }
587
588 if (forwardKEDiff < 0.0) {
589 if (forwardKEDiff + backwardKEDiff > 0.0) {
590 forwardEnergy = forwardKinetic;
591 backwardEnergy += forwardKEDiff;
592 backwardKEDiff = backwardEnergy - backwardKinetic;
593 forwardKEDiff = 0.0;
594 } else {
595 G4cout << " False return due to insufficient forward energy " << G4endl;
596 return false;
597 }
598 }
599 }
600
601 // Generate momentum for incident (current) particle, which was placed
602 // in the forward hemisphere.
603 // Set mass parameter for lambda fragmentation model.
604 // Set pt and phi values, which are changed somewhat in the iteration loop
605
606 G4double ran = -G4Log(1.0-G4UniformRand());
607 if (currentParticle.GetDefinition()->GetParticleSubType() == "pi") {
608 aspar = 0.60;
609 pt = std::sqrt( std::pow( ran/6.0, 1.7 ) );
610 } else if (currentParticle.GetDefinition()->GetParticleSubType() == "kaon") {
611 aspar = 0.50;
612 pt = std::sqrt( std::pow( ran/5.0, 1.4 ) );
613 } else {
614 aspar = 0.40;
615 pt = std::sqrt( std::pow( ran/4.0, 1.2 ) );
616 }
617
618 phi = G4UniformRand()*twopi;
619 currentParticle.SetMomentum(pt*std::cos(phi)*GeV, pt*std::sin(phi)*GeV);
620 et = pseudoParticle[0].GetTotalEnergy()/GeV;
621 dndl[0] = 0.0;
622 vecMass = currentParticle.GetMass()/GeV;
623
624 FragmentationIntegral(pt, et, aspar, vecMass);
625
626 ran = G4UniformRand()*dndl[19];
627 l = 1;
628
629 G4int loop = 0;
631 ed << " While count exceeded " << G4endl;
632 while( ( ran > dndl[l] ) && ( l < 19 ) ) { /* Loop checking, 01.09.2015, D.Wright */
633 l++;
634 loop++;
635 if (loop > 1000) {
636 G4Exception("G4RPGFragmentation::ReactionStage()", "HAD_RPG_100", JustWarning, ed);
637 break;
638 }
639 }
640
641 x = (G4double(l-1) + G4UniformRand())/19.;
642 currentParticle.SetMomentum( x*et*GeV ); // set the z-momentum
643
644 if (forwardEnergy < forwardKinetic) {
645 totalEnergy = vecMass + 0.04*std::fabs(normal());
646 G4cout << " Not enough forward energy: forwardEnergy = "
647 << forwardEnergy << " forwardKinetic = "
648 << forwardKinetic << " total energy left = "
649 << backwardKEDiff + forwardKEDiff << G4endl;
650 } else {
651 totalEnergy = vecMass + forwardEnergy - forwardKinetic;
652 forwardKinetic = forwardEnergy;
653 }
654 currentParticle.SetTotalEnergy( totalEnergy*GeV );
655 pp = std::sqrt(std::abs( totalEnergy*totalEnergy - vecMass*vecMass) )*GeV;
656 pp1 = currentParticle.GetMomentum().mag();
657
658 if (pp1 < 1.0e-6*GeV) {
659 G4ThreeVector iso = Isotropic(pp);
660 currentParticle.SetMomentum( iso.x(), iso.y(), iso.z() );
661 } else {
662 currentParticle.SetMomentum( currentParticle.GetMomentum() * (pp/pp1) );
663 }
664 pseudoParticle[4] = pseudoParticle[4] + currentParticle;
665
666 // Current particle now finished
667
668 // Begin target particle
669
670 if (backwardNucleonCount < 18) {
671 targetParticle.SetSide(-3);
672 ++backwardNucleonCount;
673
674 } else {
675 // Set pt and phi values, they are changed somewhat in the iteration loop
676 // Set mass parameter for lambda fragmentation model
677
678 vecMass = targetParticle.GetMass()/GeV;
679 ran = -G4Log(1.0-G4UniformRand());
680 aspar = 0.40;
681 pt = std::max( 0.001, std::sqrt( std::pow( ran/4.0, 1.2 ) ) );
682 phi = G4UniformRand()*twopi;
683 targetParticle.SetMomentum(pt*std::cos(phi)*GeV, pt*std::sin(phi)*GeV);
684 et = pseudoParticle[1].GetTotalEnergy()/GeV;
685 outerCounter = 0;
686 innerCounter = 0;
687 G4bool marginalEnergy = true;
688 dndl[0] = 0.0;
689 G4double xxx = 0.999;
690 if( extraNucleonCount < 20 ) xxx = 0.95+0.05*extraNucleonCount/20.0;
691 G4ThreeVector momentum;
692
693 while (++outerCounter < 4) { /* Loop checking, 01.09.2015, D.Wright */
694 FragmentationIntegral(pt, et, aspar, vecMass);
695
696 for (innerCounter = 0; innerCounter < 6; innerCounter++) {
697 ran = G4UniformRand()*dndl[19];
698 l = 1;
699
700 G4int loopa = 0;
702 eda << " While count exceeded " << G4endl;
703 while( ( ran > dndl[l] ) && ( l < 19 ) ) { /* Loop checking, 01.09.2015, D.Wright */
704 l++;
705 loopa++;
706 if (loopa > 1000) {
707 G4Exception("G4RPGFragmentation::ReactionStage()", "HAD_RPG_100", JustWarning, eda);
708 break;
709 }
710 }
711
712 x = -(G4double(l-1) + G4UniformRand())/19.;
713 targetParticle.SetMomentum( x*et*GeV ); // set the z-momentum
714 totalEnergy = std::sqrt(x*et*x*et + pt*pt + vecMass*vecMass);
715 targetParticle.SetTotalEnergy( totalEnergy*GeV );
716
717 if ((backwardKinetic+totalEnergy-vecMass) < xxx*backwardEnergy) {
718 pseudoParticle[5] = pseudoParticle[5] + targetParticle;
719 backwardKinetic += totalEnergy - vecMass;
720 outerCounter = 3; // leave outer loop
721 marginalEnergy = false;
722 break; // leave inner loop
723 }
724 momentum = targetParticle.GetMomentum();
725 targetParticle.SetMomentum(momentum.x() * 0.9, momentum.y() * 0.9);
726 pt *= 0.9;
727 dndl[19] *= 0.9;
728 }
729 }
730
731 if (marginalEnergy) {
732 G4cout << " Extra backward kinetic energy = "
733 << 0.999*backwardEnergy - backwardKinetic << G4endl;
734 totalEnergy = vecMass + 0.999*backwardEnergy - backwardKinetic;
735 targetParticle.SetTotalEnergy(totalEnergy*GeV);
736 pp = std::sqrt(std::abs(totalEnergy*totalEnergy - vecMass*vecMass) )*GeV;
737 targetParticle.SetMomentum(momentum.x()/0.9, momentum.y()/0.9);
738 pp1 = targetParticle.GetMomentum().mag();
739 targetParticle.SetMomentum(targetParticle.GetMomentum() * pp/pp1 );
740 pseudoParticle[5] = pseudoParticle[5] + targetParticle;
741 backwardKinetic = 0.999*backwardEnergy;
742 }
743
744 }
745
746 if (backwardEnergy < backwardKinetic)
747 G4cout << " Backward Edif = " << backwardEnergy - backwardKinetic << G4endl;
748 if (forwardEnergy != forwardKinetic)
749 G4cout << " Forward Edif = " << forwardEnergy - forwardKinetic << G4endl;
750
751 // Target particle finished.
752 // Now produce backward nucleons with a cluster model
753 // ps[2] = CM frame of projectile + target
754 // ps[3] = sum of projectile + nucleon cluster in lab frame
755 // ps[6] = proj + cluster 4-vector boosted into CM frame of proj + targ
756 // with secondaries, current and target particles subtracted
757 // = total 4-momentum of backward nucleon cluster
758
759 pseudoParticle[6].Lorentz( pseudoParticle[3], pseudoParticle[2] );
760 pseudoParticle[6] = pseudoParticle[6] - pseudoParticle[4];
761 pseudoParticle[6] = pseudoParticle[6] - pseudoParticle[5];
762
763 if (backwardNucleonCount == 1) {
764 // Target particle is the only backward nucleon. Give it the remainder
765 // of the backward kinetic energy.
766
767 G4double ekin =
768 std::min(backwardEnergy-backwardKinetic, centerofmassEnergy/2.0-protonMass/GeV);
769
770 if( ekin < 0.04 )ekin = 0.04 * std::fabs( normal() );
771 vecMass = targetParticle.GetMass()/GeV;
772 totalEnergy = ekin + vecMass;
773 targetParticle.SetTotalEnergy( totalEnergy*GeV );
774 pp = std::sqrt(std::abs(totalEnergy*totalEnergy - vecMass*vecMass) )*GeV;
775 pp1 = pseudoParticle[6].GetMomentum().mag();
776 if (pp1 < 1.0e-6*GeV) {
777 G4ThreeVector iso = Isotropic(pp);
778 targetParticle.SetMomentum( iso.x(), iso.y(), iso.z() );
779 } else {
780 targetParticle.SetMomentum( pseudoParticle[6].GetMomentum() * (pp/pp1));
781 }
782 pseudoParticle[5] = pseudoParticle[5] + targetParticle;
783
784 } else if (backwardNucleonCount > 1) {
785 // Share remaining energy with up to 17 backward nucleons
786
787 G4int tempCount = 5;
788 if (backwardNucleonCount < 5) tempCount = backwardNucleonCount;
789 tempCount -= 2;
790
791 G4double clusterMass = 0.;
792 if (targetParticle.GetSide() == -3)
793 clusterMass = targetParticle.GetMass()/GeV;
794 for (i = 0; i < vecLen; ++i)
795 if (vec[i]->GetSide() == -3) clusterMass += vec[i]->GetMass()/GeV;
796 clusterMass += backwardEnergy - backwardKinetic;
797
798 totalEnergy = pseudoParticle[6].GetTotalEnergy()/GeV;
799 pseudoParticle[6].SetMass(clusterMass*GeV);
800
801 pp = std::sqrt(std::abs(totalEnergy*totalEnergy -
802 clusterMass*clusterMass) )*GeV;
803 pp1 = pseudoParticle[6].GetMomentum().mag();
804 if (pp1 < 1.0e-6*GeV) {
805 G4ThreeVector iso = Isotropic(pp);
806 pseudoParticle[6].SetMomentum(iso.x(), iso.y(), iso.z());
807 } else {
808 pseudoParticle[6].SetMomentum(pseudoParticle[6].GetMomentum() * (-pp/pp1));
809 }
810
811 std::vector<G4ReactionProduct*> tempList; // ptrs to backward nucleons
812 if (targetParticle.GetSide() == -3) tempList.push_back(&targetParticle);
813 for (i = 0; i < vecLen; ++i)
814 if (vec[i]->GetSide() == -3) tempList.push_back(vec[i]);
815
816 constantCrossSection = true;
817
818 if (tempList.size() > 1) {
819 G4int n_entry = 0;
820// wgt = GenerateNBodyEventT(pseudoParticle[6].GetMass(),
821// constantCrossSection, tempList);
822// DHW: wgt not used before before being overwritten; remove it here
823 GenerateNBodyEventT(pseudoParticle[6].GetMass(), constantCrossSection,
824 tempList);
825
826 if (targetParticle.GetSide() == -3) {
827 targetParticle = *tempList[0];
828 targetParticle.Lorentz(targetParticle, pseudoParticle[6]);
829 n_entry++;
830 }
831
832 for (i = 0; i < vecLen; ++i) {
833 if (vec[i]->GetSide() == -3) {
834 *vec[i] = *tempList[n_entry];
835 vec[i]->Lorentz(*vec[i], pseudoParticle[6]);
836 n_entry++;
837 }
838 }
839 }
840 } else return false;
841
842 if (vecLen == 0) return false; // all the secondaries have been eliminated
843
844 // Lorentz transformation to lab frame
845
846 currentParticle.Lorentz( currentParticle, pseudoParticle[1] );
847 targetParticle.Lorentz( targetParticle, pseudoParticle[1] );
848 for (i = 0; i < vecLen; ++i) vec[i]->Lorentz(*vec[i], pseudoParticle[1]);
849
850 // Set leading strange particle flag.
851 // leadFlag will be true if original particle and incident particle are
852 // both strange, in which case the incident particle becomes the leading
853 // particle.
854 // leadFlag will also be true if the target particle is strange, but the
855 // incident particle is not, in which case the target particle becomes the
856 // leading particle.
857
858 G4bool leadingStrangeParticleHasChanged = true;
859 if (leadFlag)
860 {
861 if (currentParticle.GetDefinition() == leadingStrangeParticle.GetDefinition())
862 leadingStrangeParticleHasChanged = false;
863 if (leadingStrangeParticleHasChanged &&
864 (targetParticle.GetDefinition() == leadingStrangeParticle.GetDefinition()) )
865 leadingStrangeParticleHasChanged = false;
866 if( leadingStrangeParticleHasChanged )
867 {
868 for( i=0; i<vecLen; i++ )
869 {
870 if( vec[i]->GetDefinition() == leadingStrangeParticle.GetDefinition() )
871 {
872 leadingStrangeParticleHasChanged = false;
873 break;
874 }
875 }
876 }
877 if( leadingStrangeParticleHasChanged )
878 {
879 G4bool leadTest =
880 (leadingStrangeParticle.GetDefinition()->GetParticleSubType() == "kaon" ||
881 leadingStrangeParticle.GetDefinition()->GetParticleSubType() == "pi");
882 G4bool targetTest =
883 (targetParticle.GetDefinition()->GetParticleSubType() == "kaon" ||
884 targetParticle.GetDefinition()->GetParticleSubType() == "pi");
885
886 // following modified by JLC 22-Oct-97
887
888 if( (leadTest&&targetTest) || !(leadTest||targetTest) ) // both true or both false
889 {
890 targetParticle.SetDefinitionAndUpdateE( leadingStrangeParticle.GetDefinition() );
891 targetHasChanged = true;
892 }
893 else
894 {
895 currentParticle.SetDefinitionAndUpdateE( leadingStrangeParticle.GetDefinition() );
896 incidentHasChanged = false;
897 }
898 }
899 } // end of if( leadFlag )
900
901 // Get number of final state nucleons and nucleons remaining in
902 // target nucleus
903
904 std::pair<G4int, G4int> finalStateNucleons =
905 GetFinalStateNucleons(originalTarget, vec, vecLen);
906
907 G4int protonsInFinalState = finalStateNucleons.first;
908 G4int neutronsInFinalState = finalStateNucleons.second;
909
910 G4int numberofFinalStateNucleons =
911 protonsInFinalState + neutronsInFinalState;
912
913 if (currentParticle.GetDefinition()->GetBaryonNumber() == 1 &&
914 targetParticle.GetDefinition()->GetBaryonNumber() == 1 &&
915 originalIncident->GetDefinition()->GetPDGMass() <
917 numberofFinalStateNucleons++;
918
919 numberofFinalStateNucleons = std::max(1, numberofFinalStateNucleons);
920
921 G4int PinNucleus = std::max(0,
922 G4int(targetNucleus.GetZ_asInt()) - protonsInFinalState);
923 G4int NinNucleus = std::max(0,
924 G4int(targetNucleus.GetA_asInt()-targetNucleus.GetZ_asInt()) - neutronsInFinalState);
925
926 pseudoParticle[3].SetMomentum( 0.0, 0.0, pOriginal*GeV );
927 pseudoParticle[3].SetMass( mOriginal*GeV );
928 pseudoParticle[3].SetTotalEnergy(
929 std::sqrt( pOriginal*pOriginal + mOriginal*mOriginal )*GeV );
930
931 const G4ParticleDefinition* aOrgDef = modifiedOriginal.GetDefinition();
932 G4int diff = 0;
933 if(aOrgDef == G4Proton::Proton() || aOrgDef == G4Neutron::Neutron() ) diff = 1;
934 if(numberofFinalStateNucleons == 1) diff = 0;
935 pseudoParticle[4].SetMomentum( 0.0, 0.0, 0.0 );
936 pseudoParticle[4].SetMass( protonMass*(numberofFinalStateNucleons-diff) );
937 pseudoParticle[4].SetTotalEnergy( protonMass*(numberofFinalStateNucleons-diff) );
938
939 G4double theoreticalKinetic =
940 pseudoParticle[3].GetTotalEnergy() + pseudoParticle[4].GetTotalEnergy() -
941 currentParticle.GetMass() - targetParticle.GetMass();
942 for (i = 0; i < vecLen; ++i) theoreticalKinetic -= vec[i]->GetMass();
943
944 G4double simulatedKinetic =
945 currentParticle.GetKineticEnergy() + targetParticle.GetKineticEnergy();
946 for (i = 0; i < vecLen; ++i)
947 simulatedKinetic += vec[i]->GetKineticEnergy();
948
949 pseudoParticle[5] = pseudoParticle[3] + pseudoParticle[4];
950 pseudoParticle[3].Lorentz( pseudoParticle[3], pseudoParticle[5] );
951 pseudoParticle[4].Lorentz( pseudoParticle[4], pseudoParticle[5] );
952
953 pseudoParticle[7].SetZero();
954 pseudoParticle[7] = pseudoParticle[7] + currentParticle;
955 pseudoParticle[7] = pseudoParticle[7] + targetParticle;
956 for (i = 0; i < vecLen; ++i)
957 pseudoParticle[7] = pseudoParticle[7] + *vec[i];
958
959 /*
960 // This code does not appear to do anything to the energy or angle spectra
961 if( vecLen <= 16 && vecLen > 0 )
962 {
963 // must create a new set of ReactionProducts here because GenerateNBody will
964 // modify the momenta for the particles, and we don't want to do this
965 //
966 G4ReactionProduct tempR[130];
967 tempR[0] = currentParticle;
968 tempR[1] = targetParticle;
969 for( i=0; i<vecLen; ++i )tempR[i+2] = *vec[i];
970 G4FastVector<G4ReactionProduct,256> tempV1;
971 tempV1.Initialize( vecLen+2 );
972 G4int tempLen1 = 0;
973 for( i=0; i<vecLen+2; ++i )tempV1.SetElement( tempLen1++, &tempR[i] );
974 constantCrossSection = true;
975
976 wgt = GenerateNBodyEvent(pseudoParticle[3].GetTotalEnergy() +
977 pseudoParticle[4].GetTotalEnergy(),
978 constantCrossSection, tempV1, tempLen1);
979 if (wgt == -1) {
980 G4double Qvalue = 0;
981 for (i = 0; i < tempLen1; i++) Qvalue += tempV1[i]->GetMass();
982 wgt = GenerateNBodyEvent(Qvalue,
983 constantCrossSection, tempV1, tempLen1);
984 }
985 if(wgt>-.5)
986 {
987 theoreticalKinetic = 0.0;
988 for( i=0; i<tempLen1; ++i )
989 {
990 pseudoParticle[6].Lorentz( *tempV1[i], pseudoParticle[4] );
991 theoreticalKinetic += pseudoParticle[6].GetKineticEnergy();
992 }
993 }
994 // DEBUGGING --> DumpFrames::DumpFrame(vec, vecLen);
995 }
996 */
997
998 //
999 // Make sure that the kinetic energies are correct
1000 //
1001
1002 if (simulatedKinetic != 0.0) {
1003 wgt = theoreticalKinetic/simulatedKinetic;
1004 theoreticalKinetic = currentParticle.GetKineticEnergy() * wgt;
1005 simulatedKinetic = theoreticalKinetic;
1006 currentParticle.SetKineticEnergy(theoreticalKinetic);
1007 pp = currentParticle.GetTotalMomentum();
1008 pp1 = currentParticle.GetMomentum().mag();
1009 if (pp1 < 1.0e-6*GeV) {
1010 G4ThreeVector iso = Isotropic(pp);
1011 currentParticle.SetMomentum( iso.x(), iso.y(), iso.z() );
1012 } else {
1013 currentParticle.SetMomentum(currentParticle.GetMomentum() * (pp/pp1));
1014 }
1015
1016 theoreticalKinetic = targetParticle.GetKineticEnergy() * wgt;
1017 targetParticle.SetKineticEnergy(theoreticalKinetic);
1018 simulatedKinetic += theoreticalKinetic;
1019 pp = targetParticle.GetTotalMomentum();
1020 pp1 = targetParticle.GetMomentum().mag();
1021
1022 if (pp1 < 1.0e-6*GeV) {
1023 G4ThreeVector iso = Isotropic(pp);
1024 targetParticle.SetMomentum(iso.x(), iso.y(), iso.z() );
1025 } else {
1026 targetParticle.SetMomentum(targetParticle.GetMomentum() * (pp/pp1) );
1027 }
1028
1029 for (i = 0; i < vecLen; ++i ) {
1030 theoreticalKinetic = vec[i]->GetKineticEnergy() * wgt;
1031 simulatedKinetic += theoreticalKinetic;
1032 vec[i]->SetKineticEnergy(theoreticalKinetic);
1033 pp = vec[i]->GetTotalMomentum();
1034 pp1 = vec[i]->GetMomentum().mag();
1035 if( pp1 < 1.0e-6*GeV ) {
1036 G4ThreeVector iso = Isotropic(pp);
1037 vec[i]->SetMomentum(iso.x(), iso.y(), iso.z() );
1038 } else {
1039 vec[i]->SetMomentum(vec[i]->GetMomentum() * (pp/pp1) );
1040 }
1041 }
1042 }
1043
1044 // Rotate(numberofFinalStateNucleons, pseudoParticle[3].GetMomentum(),
1045 // modifiedOriginal, originalIncident, targetNucleus,
1046 // currentParticle, targetParticle, vec, vecLen );
1047
1048 // Add black track particles
1049 // the total number of particles produced is restricted to 198
1050 // this may have influence on very high energies
1051
1052 if( atomicWeight >= 1.5 )
1053 {
1054 // npnb is number of proton/neutron black track particles
1055 // ndta is the number of deuterons, tritons, and alphas produced
1056 // epnb is the kinetic energy available for proton/neutron black track
1057 // particles
1058 // edta is the kinetic energy available for deuteron/triton/alpha particles
1059
1060 G4int npnb = 0;
1061 G4int ndta = 0;
1062
1063 G4double epnb, edta;
1064 if (veryForward) {
1065 epnb = targetNucleus.GetAnnihilationPNBlackTrackEnergy();
1066 edta = targetNucleus.GetAnnihilationDTABlackTrackEnergy();
1067 } else {
1068 epnb = targetNucleus.GetPNBlackTrackEnergy();
1069 edta = targetNucleus.GetDTABlackTrackEnergy();
1070 }
1071 /*
1072 G4ReactionProduct* fudge = new G4ReactionProduct();
1073 fudge->SetDefinition( aProton );
1074 G4double TT = epnb + edta;
1075 G4double MM = fudge->GetMass()/GeV;
1076 fudge->SetTotalEnergy(MM*GeV + TT*GeV);
1077 G4double pzz = std::sqrt(TT*(TT + 2.*MM));
1078 fudge->SetMomentum(0.0, 0.0, pzz*GeV);
1079 vec.SetElement(vecLen++, fudge);
1080 // G4cout << " Fudge = " << vec[vecLen-1]->GetKineticEnergy()/GeV
1081 << G4endl;
1082 */
1083
1084 const G4double pnCutOff = 0.001;
1085 const G4double dtaCutOff = 0.001;
1086 // const G4double kineticMinimum = 1.e-6;
1087 // const G4double kineticFactor = -0.010;
1088 // G4double sprob = 0.0; // sprob = probability of self-absorption in
1089 // heavy molecules
1090 // Not currently used (DHW 9 June 2008) const G4double ekIncident = originalIncident->GetKineticEnergy()/GeV;
1091 // if (ekIncident >= 5.0) sprob = std::min(1.0, 0.6*std::log(ekIncident-4.0));
1092 if (epnb > pnCutOff)
1093 {
1094 npnb = G4Poisson((1.5+1.25*numberofFinalStateNucleons)*epnb/(epnb+edta));
1095 if (numberofFinalStateNucleons + npnb > atomicWeight)
1096 npnb = G4int(atomicWeight+0.00001 - numberofFinalStateNucleons);
1097 npnb = std::min( npnb, 127-vecLen );
1098 }
1099 if( edta >= dtaCutOff )
1100 {
1101 ndta = G4Poisson((1.5+1.25*numberofFinalStateNucleons)*edta/(epnb+edta));
1102 ndta = std::min( ndta, 127-vecLen );
1103 }
1104 if (npnb == 0 && ndta == 0) npnb = 1;
1105
1106 AddBlackTrackParticles(epnb, npnb, edta, ndta, modifiedOriginal,
1107 PinNucleus, NinNucleus, targetNucleus,
1108 vec, vecLen);
1109 }
1110
1111 // if( centerofmassEnergy <= (4.0+G4UniformRand()) )
1112 // MomentumCheck( modifiedOriginal, currentParticle, targetParticle,
1113 // vec, vecLen );
1114 //
1115 // calculate time delay for nuclear reactions
1116 //
1117
1118 if( (atomicWeight >= 1.5) && (atomicWeight <= 230.0) && (ekOriginal <= 0.2) )
1119 currentParticle.SetTOF(
1120 1.0-500.0*G4Exp(-ekOriginal/0.04)*G4Log(G4UniformRand()) );
1121 else
1122 currentParticle.SetTOF( 1.0 );
1123 return true;
1124
1125}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
G4double G4Log(G4double x)
Definition: G4Log.hh:226
G4long G4Poisson(G4double mean)
Definition: G4Poisson.hh:50
bool G4bool
Definition: G4Types.hh:86
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
double z() const
double x() const
double y() const
double mag() const
void SetElement(G4int anIndex, Type *anElement)
Definition: G4FastVector.hh:72
const G4ParticleDefinition * GetDefinition() const
static G4Lambda * Lambda()
Definition: G4Lambda.cc:107
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
G4double GetAnnihilationPNBlackTrackEnergy() const
Definition: G4Nucleus.hh:156
G4double GetPNBlackTrackEnergy() const
Definition: G4Nucleus.hh:150
G4double GetAnnihilationDTABlackTrackEnergy() const
Definition: G4Nucleus.hh:159
G4double GetDTABlackTrackEnergy() const
Definition: G4Nucleus.hh:153
const G4String & GetParticleSubType() const
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:97
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:97
static G4PionZero * PionZero()
Definition: G4PionZero.cc:107
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
G4double A13(G4double A) const
Definition: G4Pow.cc:120
static G4Proton * Proton()
Definition: G4Proton.cc:92
void FragmentationIntegral(G4double, G4double, G4double, G4double)
G4double GenerateNBodyEventT(const G4double totalEnergy, const G4bool constantCrossSection, std::vector< G4ReactionProduct * > &list)
std::pair< G4int, G4int > GetFinalStateNucleons(const G4DynamicParticle *originalTarget, const G4FastVector< G4ReactionProduct, 256 > &vec, const G4int &vecLen)
void AddBlackTrackParticles(const G4double, const G4int, const G4double, const G4int, const G4ReactionProduct &, G4int, G4int, const G4Nucleus &, G4FastVector< G4ReactionProduct, 256 > &, G4int &)
G4ThreeVector Isotropic(const G4double &)
G4double normal()
void SetMomentum(const G4double x, const G4double y, const G4double z)
void SetTotalEnergy(const G4double en)
G4double GetTotalMomentum() const
G4double GetKineticEnergy() const
const G4ParticleDefinition * GetDefinition() const
G4double GetTotalEnergy() const
void SetNewlyAdded(const G4bool f)
G4ThreeVector GetMomentum() const
void SetSide(const G4int sid)
void Lorentz(const G4ReactionProduct &p1, const G4ReactionProduct &p2)
void SetDefinitionAndUpdateE(const G4ParticleDefinition *aParticleDefinition)
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
void SetKineticEnergy(const G4double en)
G4double GetMass() const
void SetMass(const G4double mas)

Referenced by G4RPGInelastic::CalculateMomenta().


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