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
G4AntiProtonAnnihilationAtRest.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// G4AntiProtonAnnihilationAtRest 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"
38#include "Randomize.hh"
41
42#define MAX_SECONDARIES 100
43
44// constructor
45
46G4AntiProtonAnnihilationAtRest::G4AntiProtonAnnihilationAtRest(const G4String& processName,
47 G4ProcessType aType ) :
48 G4VRestProcess (processName, aType), // initialization
49 massPionMinus(G4PionMinus::PionMinus()->GetPDGMass()/GeV),
50 massProton(G4Proton::Proton()->GetPDGMass()/GeV),
51 massPionZero(G4PionZero::PionZero()->GetPDGMass()/GeV),
52 massAntiProton(G4AntiProton::AntiProton()->GetPDGMass()/GeV),
53 massPionPlus(G4PionPlus::PionPlus()->GetPDGMass()/GeV),
54 massGamma(G4Gamma::Gamma()->GetPDGMass()/GeV),
55 pdefGamma(G4Gamma::Gamma()),
56 pdefPionPlus(G4PionPlus::PionPlus()),
57 pdefPionZero(G4PionZero::PionZero()),
58 pdefPionMinus(G4PionMinus::PionMinus()),
59 pdefProton(G4Proton::Proton()),
60 pdefAntiProton(G4AntiProton::AntiProton()),
61 pdefNeutron(G4Neutron::Neutron()),
62 pdefDeuteron(G4Deuteron::Deuteron()),
63 pdefTriton(G4Triton::Triton()),
64 pdefAlpha(G4Alpha::Alpha())
65{
66 G4HadronicDeprecate("G4AntiProtonAnnihilationAtRest");
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 == pdefAntiProton );
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 << "G4AntiProtonAnnihilationAtRestProcess::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 AntiProtons at rest; a AntiProton can either create secondaries or
154// 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
198 if (verboseLevel>1) {
199 G4cout << "G4AntiProtonAnnihilationAtRest::AtRestDoIt is invoked " <<G4endl;
200 }
201
202 G4ParticleMomentum momentum;
203 G4float localtime;
204
206
207 GenerateSecondaries(); // Generate secondaries
208
210
211 for ( G4int isec = 0; isec < ngkine; isec++ ) {
212 G4DynamicParticle* aNewParticle = new G4DynamicParticle;
213 aNewParticle->SetDefinition( gkin[isec].GetParticleDef() );
214 aNewParticle->SetMomentum( gkin[isec].GetMomentum() * GeV );
215
216 localtime = globalTime + gkin[isec].GetTOF();
217
218 G4Track* aNewTrack = new G4Track( aNewParticle, localtime*s, position );
219 aNewTrack->SetTouchableHandle(track.GetTouchableHandle());
220 aParticleChange.AddSecondary( aNewTrack );
221
222 }
223
225
226 aParticleChange.ProposeTrackStatus(fStopAndKill); // Kill the incident AntiProton
227
228// clear InteractionLengthLeft
229
231
232 return &aParticleChange;
233
234}
235
236
237void G4AntiProtonAnnihilationAtRest::GenerateSecondaries()
238{
239 static G4int index;
240 static G4int l;
241 static G4int nopt;
242 static G4int i;
243 // DHW 15 May 2011: unused: static G4ParticleDefinition* jnd;
244
245 for (i = 1; i <= MAX_SECONDARIES; ++i) {
246 pv[i].SetZero();
247 }
248
249 ngkine = 0; // number of generated secondary particles
250 ntot = 0;
251 result.SetZero();
252 result.SetMass( massAntiProton );
253 result.SetKineticEnergyAndUpdate( 0. );
254 result.SetTOF( 0. );
255 result.SetParticleDef( pdefAntiProton );
256
257 AntiProtonAnnihilation(&nopt);
258
259 // *** CHECK WHETHER THERE ARE NEW PARTICLES GENERATED ***
260 if (ntot != 0 || result.GetParticleDef() != pdefAntiProton) {
261 // *** CURRENT PARTICLE IS NOT THE SAME AS IN THE BEGINNING OR/AND ***
262 // *** ONE OR MORE SECONDARIES HAVE BEEN GENERATED ***
263
264 // --- INITIAL PARTICLE TYPE HAS BEEN CHANGED ==> PUT NEW TYPE ON ---
265 // --- THE GEANT TEMPORARY STACK ---
266
267 // --- PUT PARTICLE ON THE STACK ---
268 gkin[0] = result;
269 gkin[0].SetTOF( result.GetTOF() * 5e-11 );
270 ngkine = 1;
271
272 // --- ALL QUANTITIES ARE TAKEN FROM THE GHEISHA STACK WHERE THE ---
273 // --- CONVENTION IS THE FOLLOWING ---
274
275 // --- ONE OR MORE SECONDARIES HAVE BEEN GENERATED ---
276 for (l = 1; l <= ntot; ++l) {
277 index = l - 1;
278 // DHW 15 May 2011: unused: jnd = eve[index].GetParticleDef();
279
280 // --- ADD PARTICLE TO THE STACK IF STACK NOT YET FULL ---
281 if (ngkine < MAX_SECONDARIES) {
282 gkin[ngkine] = eve[index];
283 gkin[ngkine].SetTOF( eve[index].GetTOF() * 5e-11 );
284 ++ngkine;
285 }
286 }
287 }
288 else {
289 // --- NO SECONDARIES GENERATED AND PARTICLE IS STILL THE SAME ---
290 // --- ==> COPY EVERYTHING BACK IN THE CURRENT GEANT STACK ---
291 ngkine = 0;
292 ntot = 0;
293 globalTime += result.GetTOF() * G4float(5e-11);
294 }
295
296 // --- LIMIT THE VALUE OF NGKINE IN CASE OF OVERFLOW ---
297 ngkine = G4int(std::min(ngkine,G4int(MAX_SECONDARIES)));
298
299} // GenerateSecondaries
300
301
302void G4AntiProtonAnnihilationAtRest::Poisso(G4float xav, G4int *iran)
303{
304 static G4int i;
305 static G4float r, p1, p2, p3;
306 static G4int fivex;
307 static G4float rr, ran, rrr, ran1;
308
309 // *** GENERATION OF POISSON DISTRIBUTION ***
310 // *** NVE 16-MAR-1988 CERN GENEVA ***
311 // ORIGIN : H.FESEFELDT (27-OCT-1983)
312
313 // --- USE NORMAL DISTRIBUTION FOR <X> > 9.9 ---
314 if (xav > G4float(9.9)) {
315 // ** NORMAL DISTRIBUTION WITH SIGMA**2 = <X>
316 Normal(&ran1);
317 ran1 = xav + ran1 * std::sqrt(xav);
318 *iran = G4int(ran1);
319 if (*iran < 0) {
320 *iran = 0;
321 }
322 }
323 else {
324 fivex = G4int(xav * G4float(5.));
325 *iran = 0;
326 if (fivex > 0) {
327 r = std::exp(-G4double(xav));
328 ran1 = G4UniformRand();
329 if (ran1 > r) {
330 rr = r;
331 for (i = 1; i <= fivex; ++i) {
332 ++(*iran);
333 if (i <= 5) {
334 rrr = std::pow(xav, G4float(i)) / NFac(i);
335 }
336 // ** STIRLING' S FORMULA FOR LARGE NUMBERS
337 if (i > 5) {
338 rrr = std::exp(i * std::log(xav) -
339 (i + G4float(.5)) * std::log(i * G4float(1.)) +
340 i - G4float(.9189385));
341 }
342 rr += r * rrr;
343 if (ran1 <= rr) {
344 break;
345 }
346 }
347 }
348 }
349 else {
350 // ** FOR VERY SMALL XAV TRY IRAN=1,2,3
351 p1 = xav * std::exp(-G4double(xav));
352 p2 = xav * p1 / G4float(2.);
353 p3 = xav * p2 / G4float(3.);
354 ran = G4UniformRand();
355 if (ran >= p3) {
356 if (ran >= p2) {
357 if (ran >= p1) {
358 *iran = 0;
359 }
360 else {
361 *iran = 1;
362 }
363 }
364 else {
365 *iran = 2;
366 }
367 }
368 else {
369 *iran = 3;
370 }
371 }
372 }
373
374} // Poisso
375
376
377G4int G4AntiProtonAnnihilationAtRest::NFac(G4int n)
378{
379 G4int ret_val;
380
381 static G4int i, j;
382
383 // *** NVE 16-MAR-1988 CERN GENEVA ***
384 // ORIGIN : H.FESEFELDT (27-OCT-1983)
385
386 ret_val = 1;
387 j = n;
388 if (j > 1) {
389 if (j > 10) {
390 j = 10;
391 }
392 for (i = 2; i <= j; ++i) {
393 ret_val *= i;
394 }
395 }
396 return ret_val;
397
398} // NFac
399
400
401void G4AntiProtonAnnihilationAtRest::Normal(G4float *ran)
402{
403 static G4int i;
404
405 // *** NVE 14-APR-1988 CERN GENEVA ***
406 // ORIGIN : H.FESEFELDT (27-OCT-1983)
407
408 *ran = G4float(-6.);
409 for (i = 1; i <= 12; ++i) {
410 *ran += G4UniformRand();
411 }
412
413} // Normal
414
415
416void G4AntiProtonAnnihilationAtRest::AntiProtonAnnihilation(G4int *nopt)
417{
418 static G4float brr[3] = { G4float(.125),G4float(.25),G4float(.5) };
419
420 G4float r__1;
421
422 static G4int i, ii, kk;
423 static G4int nt;
424 static G4float cfa, eka;
425 static G4int ika, nbl;
426 static G4float ran, pcm;
427 static G4int isw;
428 static G4float tex;
429 static G4ParticleDefinition* ipa1;
430 static G4float ran1, ran2, ekin, tkin;
431 static G4float targ;
432 static G4ParticleDefinition* inve;
433 static G4float ekin1, ekin2, black;
434 static G4float pnrat, rmnve1, rmnve2;
435 static G4float ek, en;
436
437 // *** ANTI PROTON ANNIHILATION AT REST ***
438 // *** NVE 04-MAR-1988 CERN GENEVA ***
439 // ORIGIN : H.FESEFELDT (09-JULY-1987)
440
441 // NOPT=0 NO ANNIHILATION
442 // NOPT=1 ANNIH.IN PI+ PI-
443 // NOPT=2 ANNIH.IN PI0 PI0
444 // NOPT=3 ANNIH.IN PI- PI0
445 // NOPT=4 ANNIH.IN GAMMA GAMMA
446
447 pv[1].SetZero();
448 pv[1].SetMass( massAntiProton );
449 pv[1].SetKineticEnergyAndUpdate( 0. );
450 pv[1].SetTOF( result.GetTOF() );
451 pv[1].SetParticleDef( result.GetParticleDef() );
452 isw = 1;
453 ran = G4UniformRand();
454 if (ran > brr[0]) {
455 isw = 2;
456 }
457 if (ran > brr[1]) {
458 isw = 3;
459 }
460 if (ran > brr[2]) {
461 isw = 4;
462 }
463 *nopt = isw;
464 // **
465 // ** EVAPORATION
466 // **
467 if (isw == 1) {
468 rmnve1 = massPionPlus;
469 rmnve2 = massPionMinus;
470 }
471 else if (isw == 2) {
472 rmnve1 = massPionZero;
473 rmnve2 = massPionZero;
474 }
475 else if (isw == 3) {
476 rmnve1 = massPionMinus;
477 rmnve2 = massPionZero;
478 }
479 else if (isw == 4) {
480 rmnve1 = massGamma;
481 rmnve2 = massGamma;
482 }
483 ek = massProton + massAntiProton - rmnve1 - rmnve2;
484 tkin = ExNu(ek);
485 ek -= tkin;
486 if (ek < G4float(1e-4)) {
487 ek = G4float(1e-4);
488 }
489 ek *= G4float(.5);
490 en = ek + (rmnve1 + rmnve2) * G4float(.5);
491 r__1 = en * en - rmnve1 * rmnve2;
492 pcm = r__1 > 0 ? std::sqrt(r__1) : 0;
493 pv[2].SetZero();
494 pv[2].SetMass( rmnve1 );
495 pv[3].SetZero();
496 pv[3].SetMass( rmnve2 );
497 if (isw > 3) {
498 pv[2].SetMass( 0. );
499 pv[3].SetMass( 0. );
500 }
501 pv[2].SetEnergyAndUpdate( std::sqrt(pv[2].GetMass()*pv[2].GetMass()+pcm*pcm) );
502 pv[2].SetTOF( result.GetTOF() );
503 pv[3].SetEnergy( std::sqrt(pv[3].GetMass()*pv[3].GetMass()+pcm*pcm) );
504 pv[3].SetMomentumAndUpdate( -pv[2].GetMomentum().x(), -pv[2].GetMomentum().y(), -pv[2].GetMomentum().z() );
505 pv[3].SetTOF( result.GetTOF() );
506 switch ((int)isw) {
507 case 1:
508 pv[2].SetParticleDef( pdefPionPlus );
509 pv[3].SetParticleDef( pdefPionMinus );
510 break;
511 case 2:
512 pv[2].SetParticleDef( pdefPionZero );
513 pv[3].SetParticleDef( pdefPionZero );
514 break;
515 case 3:
516 pv[2].SetParticleDef( pdefPionMinus );
517 pv[3].SetParticleDef( pdefPionZero );
518 break;
519 case 4:
520 pv[2].SetParticleDef( pdefGamma );
521 pv[3].SetParticleDef( pdefGamma );
522 break;
523 default:
524 break;
525 }
526 nt = 3;
527 if (targetAtomicMass >= G4float(1.5)) {
528 cfa = (targetAtomicMass - G4float(1.)) /
529 G4float(120.) * G4float(.025) *
530 std::exp(-G4double(targetAtomicMass - G4float(1.)) / G4float(120.));
531 targ = G4float(1.);
532 tex = evapEnergy1;
533 if (tex >= G4float(.001)) {
534 black = (targ * G4float(1.25) +
535 G4float(1.5)) * evapEnergy1 / (evapEnergy1 + evapEnergy3);
536 Poisso(black, &nbl);
537 if (G4float(G4int(targ) + nbl) > targetAtomicMass) {
538 nbl = G4int(targetAtomicMass - targ);
539 }
540 if (nt + nbl > (MAX_SECONDARIES - 2)) {
541 nbl = (MAX_SECONDARIES - 2) - nt;
542 }
543 if (nbl > 0) {
544 ekin = tex / nbl;
545 ekin2 = G4float(0.);
546 for (i = 1; i <= nbl; ++i) {
547 if (nt == (MAX_SECONDARIES - 2)) {
548 continue;
549 }
550 if (ekin2 > tex) {
551 break;
552 }
553 ran1 = G4UniformRand();
554 Normal(&ran2);
555 ekin1 = -G4double(ekin) * std::log(ran1) -
556 cfa * (ran2 * G4float(.5) + G4float(1.));
557 if (ekin1 < G4float(0.)) {
558 ekin1 = std::log(ran1) * G4float(-.01);
559 }
560 ekin1 *= G4float(1.);
561 ekin2 += ekin1;
562 if (ekin2 > tex) {
563 ekin1 = tex - (ekin2 - ekin1);
564 }
565 if (ekin1 < G4float(0.)) {
566 ekin1 = G4float(.001);
567 }
568 ipa1 = pdefNeutron;
569 pnrat = G4float(1.) - targetCharge / targetAtomicMass;
570 if (G4UniformRand() > pnrat) {
571 ipa1 = pdefProton;
572 }
573 ++nt;
574 pv[nt].SetZero();
575 pv[nt].SetMass( ipa1->GetPDGMass()/GeV );
576 pv[nt].SetKineticEnergyAndUpdate( ekin1 );
577 pv[nt].SetTOF( result.GetTOF() );
578 pv[nt].SetParticleDef( ipa1 );
579 }
580 if (targetAtomicMass >= G4float(230.) && ek <= G4float(2.)) {
581 ii = nt + 1;
582 kk = 0;
583 eka = ek;
584 if (eka > G4float(1.)) {
585 eka *= eka;
586 }
587 if (eka < G4float(.1)) {
588 eka = G4float(.1);
589 }
590 ika = G4int(G4float(3.6) / eka);
591 for (i = 1; i <= nt; ++i) {
592 --ii;
593 if (pv[ii].GetParticleDef() != pdefProton) {
594 continue;
595 }
596 ipa1 = pdefNeutron;
597 pv[ii].SetMass( ipa1->GetPDGMass()/GeV );
598 pv[ii].SetParticleDef( ipa1 );
599 ++kk;
600 if (kk > ika) {
601 break;
602 }
603 }
604 }
605 }
606 }
607 // **
608 // ** THEN ALSO DEUTERONS, TRITONS AND ALPHAS
609 // **
610 tex = evapEnergy3;
611 if (tex >= G4float(.001)) {
612 black = (targ * G4float(1.25) + G4float(1.5)) * evapEnergy3 /
613 (evapEnergy1 + evapEnergy3);
614 Poisso(black, &nbl);
615 if (nt + nbl > (MAX_SECONDARIES - 2)) {
616 nbl = (MAX_SECONDARIES - 2) - nt;
617 }
618 if (nbl > 0) {
619 ekin = tex / nbl;
620 ekin2 = G4float(0.);
621 for (i = 1; i <= nbl; ++i) {
622 if (nt == (MAX_SECONDARIES - 2)) {
623 continue;
624 }
625 if (ekin2 > tex) {
626 break;
627 }
628 ran1 = G4UniformRand();
629 Normal(&ran2);
630 ekin1 = -G4double(ekin) * std::log(ran1) -
631 cfa * (ran2 * G4float(.5) + G4float(1.));
632 if (ekin1 < G4float(0.)) {
633 ekin1 = std::log(ran1) * G4float(-.01);
634 }
635 ekin1 *= G4float(1.);
636 ekin2 += ekin1;
637 if (ekin2 > tex) {
638 ekin1 = tex - (ekin2 - ekin1);
639 }
640 if (ekin1 < G4float(0.)) {
641 ekin1 = G4float(.001);
642 }
643 ran = G4UniformRand();
644 inve = pdefDeuteron;
645 if (ran > G4float(.6)) {
646 inve = pdefTriton;
647 }
648 if (ran > G4float(.9)) {
649 inve = pdefAlpha;
650 }
651 ++nt;
652 pv[nt].SetZero();
653 pv[nt].SetMass( inve->GetPDGMass()/GeV );
654 pv[nt].SetKineticEnergyAndUpdate( ekin1 );
655 pv[nt].SetTOF( result.GetTOF() );
656 pv[nt].SetParticleDef( inve );
657 }
658 }
659 }
660 }
661 result = pv[2];
662 if (nt == 2) {
663 return;
664 }
665 for (i = 3; i <= nt; ++i) {
666 if (ntot >= MAX_SECONDARIES) {
667 return;
668 }
669 eve[ntot++] = pv[i];
670 }
671
672} // AntiProtonAnnihilation
673
674
675G4double G4AntiProtonAnnihilationAtRest::ExNu(G4float ek1)
676{
677 G4float ret_val, r__1;
678
679 static G4float cfa, gfa, ran1, ran2, ekin1, atno3;
680 static G4int magic;
681 static G4float fpdiv;
682
683 // *** NUCLEAR EVAPORATION AS FUNCTION OF ATOMIC NUMBER ATNO ***
684 // *** AND KINETIC ENERGY EKIN OF PRIMARY PARTICLE ***
685 // *** NVE 04-MAR-1988 CERN GENEVA ***
686 // ORIGIN : H.FESEFELDT (10-DEC-1986)
687
688 ret_val = G4float(0.);
689 if (targetAtomicMass >= G4float(1.5)) {
690 magic = 0;
691 if (G4int(targetCharge + G4float(.1)) == 82) {
692 magic = 1;
693 }
694 ekin1 = ek1;
695 if (ekin1 < G4float(.1)) {
696 ekin1 = G4float(.1);
697 }
698 if (ekin1 > G4float(4.)) {
699 ekin1 = G4float(4.);
700 }
701 // ** 0.35 VALUE AT 1 GEV
702 // ** 0.05 VALUE AT 0.1 GEV
703 cfa = G4float(.13043478260869565);
704 cfa = cfa * std::log(ekin1) + G4float(.35);
705 if (cfa < G4float(.15)) {
706 cfa = G4float(.15);
707 }
708 ret_val = cfa * G4float(7.716) * std::exp(-G4double(cfa));
709 atno3 = targetAtomicMass;
710 if (atno3 > G4float(120.)) {
711 atno3 = G4float(120.);
712 }
713 cfa = (atno3 - G4float(1.)) /
714 G4float(120.) * std::exp(-G4double(atno3 - G4float(1.)) / G4float(120.));
715 ret_val *= cfa;
716 r__1 = ekin1;
717 fpdiv = G4float(1.) - r__1 * r__1 * G4float(.25);
718 if (fpdiv < G4float(.5)) {
719 fpdiv = G4float(.5);
720 }
721 gfa = (targetAtomicMass - G4float(1.)) /
722 G4float(70.) * G4float(2.) *
723 std::exp(-G4double(targetAtomicMass - G4float(1.)) / G4float(70.));
724 evapEnergy1 = ret_val * fpdiv;
725 evapEnergy3 = ret_val - evapEnergy1;
726 Normal(&ran1);
727 Normal(&ran2);
728 if (magic == 1) {
729 ran1 = G4float(0.);
730 ran2 = G4float(0.);
731 }
732 evapEnergy1 *= ran1 * gfa + G4float(1.);
733 if (evapEnergy1 < G4float(0.)) {
734 evapEnergy1 = G4float(0.);
735 }
736 evapEnergy3 *= ran2 * gfa + G4float(1.);
737 if (evapEnergy3 < G4float(0.)) {
738 evapEnergy3 = G4float(0.);
739 }
740 while ((ret_val = evapEnergy1 + evapEnergy3) >= ek1) {
741 evapEnergy1 *= G4float(1.) - G4UniformRand() * G4float(.5);
742 evapEnergy3 *= G4float(1.) - G4UniformRand() * G4float(.5);
743 }
744 }
745 return ret_val;
746
747} // 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
G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)
void PreparePhysicsTable(const G4ParticleDefinition &)
G4bool IsApplicable(const G4ParticleDefinition &)
void BuildPhysicsTable(const G4ParticleDefinition &)
G4double GetMeanLifeTime(const G4Track &, G4ForceCondition *)
G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
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