Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4AntiNeutronAnnihilationAtRest.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// G4AntiNeutronAnnihilationAtRest physics process
27// Larry Felawka (TRIUMF), April 1998
28//---------------------------------------------------------------------
29
30#include <string.h>
31#include <cmath>
32#include <stdio.h>
33
35#include "G4SystemOfUnits.hh"
36#include "G4DynamicParticle.hh"
37#include "G4ParticleTypes.hh"
40#include "Randomize.hh"
41
42#define MAX_SECONDARIES 100
43
44// constructor
45
46G4AntiNeutronAnnihilationAtRest::G4AntiNeutronAnnihilationAtRest(const G4String& processName,
47 G4ProcessType aType) :
48 G4VRestProcess (processName, aType), // initialization
49 massPionMinus(G4PionMinus::PionMinus()->GetPDGMass()/GeV),
50 massPionZero(G4PionZero::PionZero()->GetPDGMass()/GeV),
51 massPionPlus(G4PionPlus::PionPlus()->GetPDGMass()/GeV),
52 massGamma(G4Gamma::Gamma()->GetPDGMass()/GeV),
53 massAntiNeutron(G4AntiNeutron::AntiNeutron()->GetPDGMass()/GeV),
54 massNeutron(G4Neutron::Neutron()->GetPDGMass()/GeV),
55 pdefGamma(G4Gamma::Gamma()),
56 pdefPionPlus(G4PionPlus::PionPlus()),
57 pdefPionZero(G4PionZero::PionZero()),
58 pdefPionMinus(G4PionMinus::PionMinus()),
59 pdefProton(G4Proton::Proton()),
60 pdefNeutron(G4Neutron::Neutron()),
61 pdefAntiNeutron(G4AntiNeutron::AntiNeutron()),
62 pdefDeuteron(G4Deuteron::Deuteron()),
63 pdefTriton(G4Triton::Triton()),
64 pdefAlpha(G4Alpha::Alpha())
65{
66 G4HadronicDeprecate("G4AntiNeutronAnnihilationAtRest");
67 if (verboseLevel>0) {
68 G4cout << GetProcessName() << " is created "<< G4endl;
69 }
74
76}
77
78// destructor
79
81{
83 delete [] pv;
84 delete [] eve;
85 delete [] gkin;
86}
87
89{
91}
92
94{
96}
97
98// methods.............................................................................
99
101 const G4ParticleDefinition& particle
102 )
103{
104 return ( &particle == pdefAntiNeutron );
105
106}
107
108// Warning - this method may be optimized away if made "inline"
110{
111 return ( ngkine );
112
113}
114
115// Warning - this method may be optimized away if made "inline"
117{
118 return ( &gkin[0] );
119
120}
121
123 const G4Track& track,
125 )
126{
127 // beggining of tracking
129
130 // condition is set to "Not Forced"
132
133 // get mean life time
135
136 if ((currentInteractionLength <0.0) || (verboseLevel>2)){
137 G4cout << "G4AntiNeutronAnnihilationAtRestProcess::AtRestGetPhysicalInteractionLength ";
138 G4cout << "[ " << GetProcessName() << "]" <<G4endl;
139 track.GetDynamicParticle()->DumpInfo();
140 G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl;
141 G4cout << "MeanLifeTime = " << currentInteractionLength/ns << "[ns]" <<G4endl;
142 }
143
145
146}
147
149 const G4Track& track,
150 const G4Step&
151 )
152//
153// Handles AntiNeutrons at rest; an AntiNeutron can either create secondaries
154// or do nothing (in which case it should be sent back to decay-handling
155// section
156//
157{
158
159// Initialize ParticleChange
160// all members of G4VParticleChange are set to equal to
161// corresponding member in G4Track
162
164
165// Store some global quantities that depend on current material and particle
166
167 globalTime = track.GetGlobalTime()/s;
168 G4Material * aMaterial = track.GetMaterial();
169 const G4int numberOfElements = aMaterial->GetNumberOfElements();
170 const G4ElementVector* theElementVector = aMaterial->GetElementVector();
171
172 const G4double* theAtomicNumberDensity = aMaterial->GetAtomicNumDensityVector();
173 G4double normalization = 0;
174 for ( G4int i1=0; i1 < numberOfElements; i1++ )
175 {
176 normalization += theAtomicNumberDensity[i1] ; // change when nucleon specific
177 // probabilities are included.
178 }
179 G4double runningSum= 0.;
180 G4double random = G4UniformRand()*normalization;
181 for ( G4int i2=0; i2 < numberOfElements; i2++ )
182 {
183 runningSum += theAtomicNumberDensity[i2]; // change when nucleon specific
184 // probabilities are included.
185 if (random<=runningSum)
186 {
187 targetCharge = G4double( ((*theElementVector)[i2])->GetZ());
188 targetAtomicMass = (*theElementVector)[i2]->GetN();
189 }
190 }
191 if (random>runningSum)
192 {
193 targetCharge = G4double( ((*theElementVector)[numberOfElements-1])->GetZ());
194 targetAtomicMass = (*theElementVector)[numberOfElements-1]->GetN();
195 }
196
197 if (verboseLevel>1) {
198 G4cout << "G4AntiNeutronAnnihilationAtRest::AtRestDoIt is invoked " <<G4endl;
199 }
200
201 G4ParticleMomentum momentum;
202 G4float localtime;
203
205
206 GenerateSecondaries(); // Generate secondaries
207
209
210 for ( G4int isec = 0; isec < ngkine; isec++ ) {
211 G4DynamicParticle* aNewParticle = new G4DynamicParticle;
212 aNewParticle->SetDefinition( gkin[isec].GetParticleDef() );
213 aNewParticle->SetMomentum( gkin[isec].GetMomentum() * GeV );
214
215 localtime = globalTime + gkin[isec].GetTOF();
216
217 G4Track* aNewTrack = new G4Track( aNewParticle, localtime*s, position );
218 aNewTrack->SetTouchableHandle(track.GetTouchableHandle());
219 aParticleChange.AddSecondary( aNewTrack );
220
221 }
222
224
225 aParticleChange.ProposeTrackStatus(fStopAndKill); // Kill the incident AntiNeutron
226
227// clear InteractionLengthLeft
228
230
231 return &aParticleChange;
232
233}
234
235
236void G4AntiNeutronAnnihilationAtRest::GenerateSecondaries()
237{
238 static G4int index;
239 static G4int l;
240 static G4int nopt;
241 static G4int i;
242 // DHW 15 May 2011: unused: static G4ParticleDefinition* jnd;
243
244 for (i = 1; i <= MAX_SECONDARIES; ++i) {
245 pv[i].SetZero();
246 }
247
248
249 ngkine = 0; // number of generated secondary particles
250 ntot = 0;
251 result.SetZero();
252 result.SetMass( massAntiNeutron );
253 result.SetKineticEnergyAndUpdate( 0. );
254 result.SetTOF( 0. );
255 result.SetParticleDef( pdefAntiNeutron );
256
257 // *** SELECT PROCESS FOR CURRENT PARTICLE ***
258
259 AntiNeutronAnnihilation(&nopt);
260
261 // *** CHECK WHETHER THERE ARE NEW PARTICLES GENERATED ***
262 if (ntot != 0 || result.GetParticleDef() != pdefAntiNeutron) {
263 // *** CURRENT PARTICLE IS NOT THE SAME AS IN THE BEGINNING OR/AND ***
264 // *** ONE OR MORE SECONDARIES HAVE BEEN GENERATED ***
265
266 // --- INITIAL PARTICLE TYPE HAS BEEN CHANGED ==> PUT NEW TYPE ON ---
267 // --- THE GEANT TEMPORARY STACK ---
268
269 // --- PUT PARTICLE ON THE STACK ---
270 gkin[0] = result;
271 gkin[0].SetTOF( result.GetTOF() * 5e-11 );
272 ngkine = 1;
273
274 // --- ALL QUANTITIES ARE TAKEN FROM THE GHEISHA STACK WHERE THE ---
275 // --- CONVENTION IS THE FOLLOWING ---
276
277 // --- ONE OR MORE SECONDARIES HAVE BEEN GENERATED ---
278 for (l = 1; l <= ntot; ++l) {
279 index = l - 1;
280 // DHW 15 May 2011: unused: jnd = eve[index].GetParticleDef();
281
282 // --- ADD PARTICLE TO THE STACK IF STACK NOT YET FULL ---
283 if (ngkine < MAX_SECONDARIES) {
284 gkin[ngkine] = eve[index];
285 gkin[ngkine].SetTOF( eve[index].GetTOF() * 5e-11 );
286 ++ngkine;
287 }
288 }
289 }
290 else {
291 // --- NO SECONDARIES GENERATED AND PARTICLE IS STILL THE SAME ---
292 // --- ==> COPY EVERYTHING BACK IN THE CURRENT GEANT STACK ---
293 ngkine = 0;
294 ntot = 0;
295 globalTime += result.GetTOF() * G4float(5e-11);
296 }
297
298 // --- LIMIT THE VALUE OF NGKINE IN CASE OF OVERFLOW ---
299 ngkine = G4int(std::min(ngkine,G4int(MAX_SECONDARIES)));
300
301} // GenerateSecondaries
302
303
304void G4AntiNeutronAnnihilationAtRest::Poisso(G4float xav, G4int *iran)
305{
306 static G4int i;
307 static G4float r, p1, p2, p3;
308 static G4int fivex;
309 static G4float rr, ran, rrr, ran1;
310
311 // *** GENERATION OF POISSON DISTRIBUTION ***
312 // *** NVE 16-MAR-1988 CERN GENEVA ***
313 // ORIGIN : H.FESEFELDT (27-OCT-1983)
314
315 // --- USE NORMAL DISTRIBUTION FOR <X> > 9.9 ---
316 if (xav > G4float(9.9)) {
317 // ** NORMAL DISTRIBUTION WITH SIGMA**2 = <X>
318 Normal(&ran1);
319 ran1 = xav + ran1 * std::sqrt(xav);
320 *iran = G4int(ran1);
321 if (*iran < 0) {
322 *iran = 0;
323 }
324 }
325 else {
326 fivex = G4int(xav * G4float(5.));
327 *iran = 0;
328 if (fivex > 0) {
329 r = std::exp(-G4double(xav));
330 ran1 = G4UniformRand();
331 if (ran1 > r) {
332 rr = r;
333 for (i = 1; i <= fivex; ++i) {
334 ++(*iran);
335 if (i <= 5) {
336 rrr = std::pow(xav, G4float(i)) / NFac(i);
337 }
338 // ** STIRLING' S FORMULA FOR LARGE NUMBERS
339 if (i > 5) {
340 rrr = std::exp(i * std::log(xav) -
341 (i + G4float(.5)) * std::log(i * G4float(1.)) +
342 i - G4float(.9189385));
343 }
344 rr += r * rrr;
345 if (ran1 <= rr) {
346 break;
347 }
348 }
349 }
350 }
351 else {
352 // ** FOR VERY SMALL XAV TRY IRAN=1,2,3
353 p1 = xav * std::exp(-G4double(xav));
354 p2 = xav * p1 / G4float(2.);
355 p3 = xav * p2 / G4float(3.);
356 ran = G4UniformRand();
357 if (ran >= p3) {
358 if (ran >= p2) {
359 if (ran >= p1) {
360 *iran = 0;
361 }
362 else {
363 *iran = 1;
364 }
365 }
366 else {
367 *iran = 2;
368 }
369 }
370 else {
371 *iran = 3;
372 }
373 }
374 }
375
376} // Poisso
377
378
379G4int G4AntiNeutronAnnihilationAtRest::NFac(G4int n)
380{
381 G4int ret_val;
382
383 static G4int i, j;
384
385 // *** NVE 16-MAR-1988 CERN GENEVA ***
386 // ORIGIN : H.FESEFELDT (27-OCT-1983)
387
388 ret_val = 1;
389 j = n;
390 if (j > 1) {
391 if (j > 10) {
392 j = 10;
393 }
394 for (i = 2; i <= j; ++i) {
395 ret_val *= i;
396 }
397 }
398 return ret_val;
399
400} // NFac
401
402
403void G4AntiNeutronAnnihilationAtRest::Normal(G4float *ran)
404{
405 static G4int i;
406
407 // *** NVE 14-APR-1988 CERN GENEVA ***
408 // ORIGIN : H.FESEFELDT (27-OCT-1983)
409
410 *ran = G4float(-6.);
411 for (i = 1; i <= 12; ++i) {
412 *ran += G4UniformRand();
413 }
414
415} // Normal
416
417
418void G4AntiNeutronAnnihilationAtRest::AntiNeutronAnnihilation(G4int *nopt)
419{
420 static G4float brr[3] = { G4float(.125),G4float(.25),G4float(.5) };
421
422 G4float r__1;
423
424 static G4int i, ii, kk;
425 static G4int nt;
426 static G4float cfa, eka;
427 static G4int ika, nbl;
428 static G4float ran, pcm;
429 static G4int isw;
430 static G4float tex;
431 static G4ParticleDefinition* ipa1;
432 static G4float ran1, ran2, ekin, tkin;
433 static G4float targ;
434 static G4ParticleDefinition* inve;
435 static G4float ekin1, ekin2, black;
436 static G4float pnrat, rmnve1, rmnve2;
437 static G4float ek, en;
438
439 // *** ANTI NEUTRON ANNIHILATION AT REST ***
440 // *** NVE 04-MAR-1988 CERN GENEVA ***
441 // ORIGIN : H.FESEFELDT (09-JULY-1987)
442
443 // NOPT=0 NO ANNIHILATION
444 // NOPT=1 ANNIH.IN PI+ PI-
445 // NOPT=2 ANNIH.IN PI0 PI0
446 // NOPT=3 ANNIH.IN PI+ PI0
447 // NOPT=4 ANNIH.IN GAMMA GAMMA
448
449 pv[1].SetZero();
450 pv[1].SetMass( massAntiNeutron );
451 pv[1].SetKineticEnergyAndUpdate( 0. );
452 pv[1].SetTOF( result.GetTOF() );
453 pv[1].SetParticleDef( result.GetParticleDef() );
454 isw = 1;
455 ran = G4UniformRand();
456 if (ran > brr[0]) {
457 isw = 2;
458 }
459 if (ran > brr[1]) {
460 isw = 3;
461 }
462 if (ran > brr[2]) {
463 isw = 4;
464 }
465 *nopt = isw;
466 // **
467 // ** EVAPORATION
468 // **
469 rmnve1 = massPionPlus;
470 rmnve2 = massPionMinus;
471 if (isw == 2) {
472 rmnve1 = massPionZero;
473 }
474 if (isw == 2) {
475 rmnve2 = massPionZero;
476 }
477 if (isw == 3) {
478 rmnve2 = massPionZero;
479 }
480 if (isw == 4) {
481 rmnve1 = massGamma;
482 }
483 if (isw == 4) {
484 rmnve2 = massGamma;
485 }
486 ek = massNeutron + massAntiNeutron - rmnve1 - rmnve2;
487 tkin = ExNu(ek);
488 ek -= tkin;
489 if (ek < G4float(1e-4)) {
490 ek = G4float(1e-4);
491 }
492 ek /= G4float(2.);
493 en = ek + (rmnve1 + rmnve2) / G4float(2.);
494 r__1 = en * en - rmnve1 * rmnve2;
495 pcm = r__1 > 0 ? std::sqrt(r__1) : 0;
496 pv[2].SetZero();
497 pv[2].SetMass( rmnve1 );
498 pv[3].SetZero();
499 pv[3].SetMass( rmnve2 );
500 if (isw > 3) {
501 pv[2].SetMass( 0. );
502 pv[3].SetMass( 0. );
503 }
504 pv[2].SetEnergyAndUpdate( std::sqrt(pv[2].GetMass()*pv[2].GetMass()+pcm*pcm) );
505 pv[2].SetTOF( result.GetTOF() );
506 pv[3].SetEnergy( std::sqrt(pv[3].GetMass()*pv[3].GetMass()+pcm*pcm) );
507 pv[3].SetMomentumAndUpdate( -pv[2].GetMomentum().x(), -pv[2].GetMomentum().y(), -pv[2].GetMomentum().z() );
508 pv[3].SetTOF( result.GetTOF() );
509 switch ((int)isw) {
510 case 1:
511 pv[2].SetParticleDef( pdefPionPlus );
512 pv[3].SetParticleDef( pdefPionMinus );
513 break;
514 case 2:
515 pv[2].SetParticleDef( pdefPionZero );
516 pv[3].SetParticleDef( pdefPionZero );
517 break;
518 case 3:
519 pv[2].SetParticleDef( pdefPionPlus );
520 pv[3].SetParticleDef( pdefPionZero );
521 break;
522 case 4:
523 pv[2].SetParticleDef( pdefGamma );
524 pv[3].SetParticleDef( pdefGamma );
525 break;
526 default:
527 break;
528 }
529 nt = 3;
530 if (targetAtomicMass >= G4float(1.5)) {
531 cfa = (targetAtomicMass - G4float(1.)) / G4float(120.) *
532 G4float(.025) * std::exp(-G4double(targetAtomicMass - G4float(1.)) /
533 G4float(120.));
534 targ = G4float(1.);
535 tex = evapEnergy1;
536 if (tex >= G4float(.001)) {
537 black = (targ * G4float(1.25) +
538 G4float(1.5)) * evapEnergy1 / (evapEnergy1 + evapEnergy3);
539 Poisso(black, &nbl);
540 if (G4float(G4int(targ) + nbl) > targetAtomicMass) {
541 nbl = G4int(targetAtomicMass - targ);
542 }
543 if (nt + nbl > (MAX_SECONDARIES - 2)) {
544 nbl = (MAX_SECONDARIES - 2) - nt;
545 }
546 if (nbl > 0) {
547 ekin = tex / nbl;
548 ekin2 = G4float(0.);
549 for (i = 1; i <= nbl; ++i) {
550 if (nt == (MAX_SECONDARIES - 2)) {
551 continue;
552 }
553 if (ekin2 > tex) {
554 break;
555 }
556 ran1 = G4UniformRand();
557 Normal(&ran2);
558 ekin1 = -G4double(ekin) * std::log(ran1) -
559 cfa * (ran2 * G4float(.5) + G4float(1.));
560 if (ekin1 < G4float(0.)) {
561 ekin1 = std::log(ran1) * G4float(-.01);
562 }
563 ekin1 *= G4float(1.);
564 ekin2 += ekin1;
565 if (ekin2 > tex) {
566 ekin1 = tex - (ekin2 - ekin1);
567 }
568 if (ekin1 < G4float(0.)) {
569 ekin1 = G4float(.001);
570 }
571 ipa1 = pdefNeutron;
572 pnrat = G4float(1.) - targetCharge / targetAtomicMass;
573 if (G4UniformRand() > pnrat) {
574 ipa1 = pdefProton;
575 }
576 ++nt;
577 pv[nt].SetZero();
578 pv[nt].SetMass( ipa1->GetPDGMass()/GeV );
579 pv[nt].SetKineticEnergyAndUpdate( ekin1 );
580 pv[nt].SetTOF( result.GetTOF() );
581 pv[nt].SetParticleDef( ipa1 );
582 }
583 if (targetAtomicMass >= G4float(230.) && ek <= G4float(2.)) {
584 ii = nt + 1;
585 kk = 0;
586 eka = ek;
587 if (eka > G4float(1.)) {
588 eka *= eka;
589 }
590 if (eka < G4float(.1)) {
591 eka = G4float(.1);
592 }
593 ika = G4int(G4float(3.6) / eka);
594 for (i = 1; i <= nt; ++i) {
595 --ii;
596 if (pv[ii].GetParticleDef() != pdefProton) {
597 continue;
598 }
599 ipa1 = pdefNeutron;
600 pv[ii].SetMass( ipa1->GetPDGMass()/GeV );
601 pv[ii].SetParticleDef( ipa1 );
602 ++kk;
603 if (kk > ika) {
604 break;
605 }
606 }
607 }
608 }
609 }
610 // **
611 // ** THEN ALSO DEUTERONS, TRITONS AND ALPHAS
612 // **
613 tex = evapEnergy3;
614 if (tex >= G4float(.001)) {
615 black = (targ * G4float(1.25) + G4float(1.5)) * evapEnergy3 /
616 (evapEnergy1 + evapEnergy3);
617 Poisso(black, &nbl);
618 if (nt + nbl > (MAX_SECONDARIES - 2)) {
619 nbl = (MAX_SECONDARIES - 2) - nt;
620 }
621 if (nbl > 0) {
622 ekin = tex / nbl;
623 ekin2 = G4float(0.);
624 for (i = 1; i <= nbl; ++i) {
625 if (nt == (MAX_SECONDARIES - 2)) {
626 continue;
627 }
628 if (ekin2 > tex) {
629 break;
630 }
631 ran1 = G4UniformRand();
632 Normal(&ran2);
633 ekin1 = -G4double(ekin) * std::log(ran1) -
634 cfa * (ran2 * G4float(.5) + G4float(1.));
635 if (ekin1 < G4float(0.)) {
636 ekin1 = std::log(ran1) * G4float(-.01);
637 }
638 ekin1 *= G4float(1.);
639 ekin2 += ekin1;
640 if (ekin2 > tex) {
641 ekin1 = tex - (ekin2 - ekin1);
642 }
643 if (ekin1 < G4float(0.)) {
644 ekin1 = G4float(.001);
645 }
646 ran = G4UniformRand();
647 inve = pdefDeuteron;
648 if (ran > G4float(.6)) {
649 inve = pdefTriton;
650 }
651 if (ran > G4float(.9)) {
652 inve = pdefAlpha;
653 }
654 ++nt;
655 pv[nt].SetZero();
656 pv[nt].SetMass( inve->GetPDGMass()/GeV );
657 pv[nt].SetKineticEnergyAndUpdate( ekin1 );
658 pv[nt].SetTOF( result.GetTOF() );
659 pv[nt].SetParticleDef( inve );
660 }
661 }
662 }
663 }
664 result = pv[2];
665 if (nt == 2) {
666 return;
667 }
668 for (i = 3; i <= nt; ++i) {
669 if (ntot >= MAX_SECONDARIES) {
670 return;
671 }
672 eve[ntot++] = pv[i];
673 }
674
675} // AntiNeutronAnnihilation
676
677
678G4double G4AntiNeutronAnnihilationAtRest::ExNu(G4float ek1)
679{
680 G4float ret_val, r__1;
681
682 static G4float cfa, gfa, ran1, ran2, ekin1, atno3;
683 static G4int magic;
684 static G4float fpdiv;
685
686 // *** NUCLEAR EVAPORATION AS FUNCTION OF ATOMIC NUMBER ATNO ***
687 // *** AND KINETIC ENERGY EKIN OF PRIMARY PARTICLE ***
688 // *** NVE 04-MAR-1988 CERN GENEVA ***
689 // ORIGIN : H.FESEFELDT (10-DEC-1986)
690
691 ret_val = G4float(0.);
692 if (targetAtomicMass >= G4float(1.5)) {
693 magic = 0;
694 if (G4int(targetCharge + G4float(.1)) == 82) {
695 magic = 1;
696 }
697 ekin1 = ek1;
698 if (ekin1 < G4float(.1)) {
699 ekin1 = G4float(.1);
700 }
701 if (ekin1 > G4float(4.)) {
702 ekin1 = G4float(4.);
703 }
704 // ** 0.35 VALUE AT 1 GEV
705 // ** 0.05 VALUE AT 0.1 GEV
706 cfa = G4float(.13043478260869565);
707 cfa = cfa * std::log(ekin1) + G4float(.35);
708 if (cfa < G4float(.15)) {
709 cfa = G4float(.15);
710 }
711 ret_val = cfa * G4float(7.716) * std::exp(-G4double(cfa));
712 atno3 = targetAtomicMass;
713 if (atno3 > G4float(120.)) {
714 atno3 = G4float(120.);
715 }
716 cfa = (atno3 - G4float(1.)) /
717 G4float(120.) * std::exp(-G4double(atno3 - G4float(1.)) / G4float(120.));
718 ret_val *= cfa;
719 r__1 = ekin1;
720 fpdiv = G4float(1.) - r__1 * r__1 * G4float(.25);
721 if (fpdiv < G4float(.5)) {
722 fpdiv = G4float(.5);
723 }
724 gfa = (targetAtomicMass - G4float(1.)) /
725 G4float(70.) * G4float(2.) *
726 std::exp(-G4double(targetAtomicMass - G4float(1.)) / G4float(70.));
727 evapEnergy1 = ret_val * fpdiv;
728 evapEnergy3 = ret_val - evapEnergy1;
729 Normal(&ran1);
730 Normal(&ran2);
731 if (magic == 1) {
732 ran1 = G4float(0.);
733 ran2 = G4float(0.);
734 }
735 evapEnergy1 *= ran1 * gfa + G4float(1.);
736 if (evapEnergy1 < G4float(0.)) {
737 evapEnergy1 = G4float(0.);
738 }
739 evapEnergy3 *= ran2 * gfa + G4float(1.);
740 if (evapEnergy3 < G4float(0.)) {
741 evapEnergy3 = G4float(0.);
742 }
743 while ((ret_val = evapEnergy1 + evapEnergy3) >= ek1) {
744 evapEnergy1 *= G4float(1.) - G4UniformRand() * G4float(.5);
745 evapEnergy3 *= G4float(1.) - G4UniformRand() * G4float(.5);
746 }
747 }
748 return ret_val;
749
750} // ExNu
#define MAX_SECONDARIES
std::vector< G4Element * > G4ElementVector
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
@ NotForced
#define G4HadronicDeprecate(name)
@ fHadronAtRest
G4ProcessType
@ fStopAndKill
double G4double
Definition: G4Types.hh:64
float G4float
Definition: G4Types.hh:65
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:53
void PreparePhysicsTable(const G4ParticleDefinition &)
G4double GetMeanLifeTime(const G4Track &, G4ForceCondition *)
G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
void BuildPhysicsTable(const G4ParticleDefinition &)
G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)
G4bool IsApplicable(const G4ParticleDefinition &)
void DumpInfo(G4int mode=0) const
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
void SetMomentum(const G4ThreeVector &momentum)
void SetEnergyAndUpdate(G4double e)
G4ParticleDefinition * GetParticleDef()
void SetParticleDef(G4ParticleDefinition *c)
void SetMomentumAndUpdate(G4ParticleMomentum mom)
void SetKineticEnergyAndUpdate(G4double ekin)
void DeRegisterExtraProcess(G4VProcess *)
void RegisterExtraProcess(G4VProcess *)
void RegisterParticleForExtraProcess(G4VProcess *, const G4ParticleDefinition *)
static G4HadronicProcessStore * Instance()
void PrintInfo(const G4ParticleDefinition *)
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:189
size_t GetNumberOfElements() const
Definition: G4Material.hh:185
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:215
const G4String & GetName() const
Definition: G4Material.hh:177
void AddSecondary(G4Track *aSecondary)
virtual void Initialize(const G4Track &)
Definition: G4Step.hh:78
const G4ThreeVector & GetPosition() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4double GetGlobalTime() const
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
const G4TouchableHandle & GetTouchableHandle() const
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetNumberOfSecondaries(G4int totSecondaries)
G4double currentInteractionLength
Definition: G4VProcess.hh:297
virtual void ResetNumberOfInteractionLengthLeft()
Definition: G4VProcess.cc:92
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
G4int verboseLevel
Definition: G4VProcess.hh:368
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:293
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:403
const G4String & GetProcessName() const
Definition: G4VProcess.hh:379
#define ns
Definition: xmlparse.cc:597