Geant4  10.02.p03
G4AntiNeutronAnnihilationAtRest Class Reference

#include <G4AntiNeutronAnnihilationAtRest.hh>

Inheritance diagram for G4AntiNeutronAnnihilationAtRest:
Collaboration diagram for G4AntiNeutronAnnihilationAtRest:

Public Member Functions

 G4AntiNeutronAnnihilationAtRest (const G4String &processName="AntiNeutronAnnihilationAtRest", G4ProcessType aType=fHadronic)
 
 ~G4AntiNeutronAnnihilationAtRest ()
 
G4bool IsApplicable (const G4ParticleDefinition &)
 
void PreparePhysicsTable (const G4ParticleDefinition &)
 
void BuildPhysicsTable (const G4ParticleDefinition &)
 
G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
G4double GetMeanLifeTime (const G4Track &, G4ForceCondition *)
 
G4VParticleChange * AtRestDoIt (const G4Track &, const G4Step &)
 
G4int GetNumberOfSecondaries ()
 
G4GHEKinematicsVectorGetSecondaryKinematics ()
 
- Public Member Functions inherited from G4VRestProcess
 G4VRestProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VRestProcess (G4VRestProcess &)
 
virtual ~G4VRestProcess ()
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &, G4double, G4ForceCondition *)
 
virtual G4VParticleChange * PostStepDoIt (const G4Track &, const G4Step &)
 
virtual G4VParticleChange * AlongStepDoIt (const G4Track &, const G4Step &)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4int operator== (const G4VProcess &right) const
 
G4int operator!= (const G4VProcess &right) const
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual void StartTracking (G4Track *)
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Private Member Functions

G4AntiNeutronAnnihilationAtRestoperator= (const G4AntiNeutronAnnihilationAtRest &right)
 
 G4AntiNeutronAnnihilationAtRest (const G4AntiNeutronAnnihilationAtRest &)
 
void GenerateSecondaries ()
 
void Poisso (G4float, G4int *)
 
void Normal (G4float *)
 
void AntiNeutronAnnihilation (G4int *)
 
G4double ExNu (G4float)
 
G4int NFac (G4int)
 

Private Attributes

G4float globalTime
 
G4float targetAtomicMass
 
G4float targetCharge
 
G4GHEKinematicsVectorpv
 
G4GHEKinematicsVectoreve
 
G4GHEKinematicsVectorgkin
 
G4float evapEnergy1
 
G4float evapEnergy3
 
G4int ngkine
 
G4int ntot
 
G4GHEKinematicsVector result
 
G4float massPionMinus
 
G4float massPionZero
 
G4float massPionPlus
 
G4float massGamma
 
G4float massAntiNeutron
 
G4float massNeutron
 
G4ParticleDefinitionpdefGamma
 
G4ParticleDefinitionpdefPionPlus
 
G4ParticleDefinitionpdefPionZero
 
G4ParticleDefinitionpdefPionMinus
 
G4ParticleDefinitionpdefProton
 
G4ParticleDefinitionpdefNeutron
 
G4ParticleDefinitionpdefAntiNeutron
 
G4ParticleDefinitionpdefDeuteron
 
G4ParticleDefinitionpdefTriton
 
G4ParticleDefinitionpdefAlpha
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager
 
G4VParticleChange * pParticleChange
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft
 
G4double currentInteractionLength
 
G4double theInitialNumberOfInteractionLength
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType
 
G4int theProcessSubType
 
G4double thePILfactor
 
G4bool enableAtRestDoIt
 
G4bool enableAlongStepDoIt
 
G4bool enablePostStepDoIt
 
G4int verboseLevel
 

Detailed Description

Definition at line 47 of file G4AntiNeutronAnnihilationAtRest.hh.

Constructor & Destructor Documentation

◆ G4AntiNeutronAnnihilationAtRest() [1/2]

G4AntiNeutronAnnihilationAtRest::G4AntiNeutronAnnihilationAtRest ( const G4AntiNeutronAnnihilationAtRest )
private

◆ G4AntiNeutronAnnihilationAtRest() [2/2]

G4AntiNeutronAnnihilationAtRest::G4AntiNeutronAnnihilationAtRest ( const G4String processName = "AntiNeutronAnnihilationAtRest",
G4ProcessType  aType = fHadronic 
)

Definition at line 46 of file G4AntiNeutronAnnihilationAtRest.cc.

47  :
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),
54  massNeutron(G4Neutron::Neutron()->GetPDGMass()/GeV),
65 {
66  G4HadronicDeprecate("G4AntiNeutronAnnihilationAtRest");
67  if (verboseLevel>0) {
68  G4cout << GetProcessName() << " is created "<< G4endl;
69  }
74 
77  = evapEnergy3 = 0.0;
78  ngkine = ntot = 0;
79 }
G4int verboseLevel
Definition: G4VProcess.hh:368
static G4HadronicProcessStore * Instance()
#define G4HadronicDeprecate(name)
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4GLOB_DLL std::ostream G4cout
#define MAX_SECONDARIES
static G4Triton * Triton()
Definition: G4Triton.cc:95
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:98
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
static const double GeV
Definition: G4SIunits.hh:214
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
static G4PionZero * PionZero()
Definition: G4PionZero.cc:108
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:94
void RegisterExtraProcess(G4VProcess *)
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:98
#define G4endl
Definition: G4ios.hh:61
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
static G4AntiNeutron * AntiNeutron()
Here is the call graph for this function:

◆ ~G4AntiNeutronAnnihilationAtRest()

G4AntiNeutronAnnihilationAtRest::~G4AntiNeutronAnnihilationAtRest ( )

Definition at line 83 of file G4AntiNeutronAnnihilationAtRest.cc.

84 {
86  delete [] pv;
87  delete [] eve;
88  delete [] gkin;
89 }
void DeRegisterExtraProcess(G4VProcess *)
static G4HadronicProcessStore * Instance()
Here is the call graph for this function:

Member Function Documentation

◆ AntiNeutronAnnihilation()

void G4AntiNeutronAnnihilationAtRest::AntiNeutronAnnihilation ( G4int nopt)
private

Definition at line 421 of file G4AntiNeutronAnnihilationAtRest.cc.

422 {
423  static G4ThreadLocal G4float brr[3] = { G4float(.125),G4float(.25),G4float(.5) };
424 
425  G4float r__1;
426 
427  static G4ThreadLocal G4int i, ii, kk;
428  static G4ThreadLocal G4int nt;
429  static G4ThreadLocal G4float cfa, eka;
430  static G4ThreadLocal G4int ika, nbl;
431  static G4ThreadLocal G4float ran, pcm;
432  static G4ThreadLocal G4int isw;
433  static G4ThreadLocal G4float tex;
434  static G4ThreadLocal G4ParticleDefinition* ipa1;
435  static G4ThreadLocal G4float ran1, ran2, ekin, tkin;
436  static G4ThreadLocal G4float targ;
437  static G4ThreadLocal G4ParticleDefinition* inve;
438  static G4ThreadLocal G4float ekin1, ekin2, black;
439  static G4ThreadLocal G4float pnrat, rmnve1, rmnve2;
440  static G4ThreadLocal G4float ek, en;
441 
442  // *** ANTI NEUTRON ANNIHILATION AT REST ***
443  // *** NVE 04-MAR-1988 CERN GENEVA ***
444  // ORIGIN : H.FESEFELDT (09-JULY-1987)
445 
446  // NOPT=0 NO ANNIHILATION
447  // NOPT=1 ANNIH.IN PI+ PI-
448  // NOPT=2 ANNIH.IN PI0 PI0
449  // NOPT=3 ANNIH.IN PI+ PI0
450  // NOPT=4 ANNIH.IN GAMMA GAMMA
451 
452  pv[1].SetZero();
453  pv[1].SetMass( massAntiNeutron );
454  pv[1].SetKineticEnergyAndUpdate( 0. );
455  pv[1].SetTOF( result.GetTOF() );
457  isw = 1;
458  ran = G4UniformRand();
459  if (ran > brr[0]) {
460  isw = 2;
461  }
462  if (ran > brr[1]) {
463  isw = 3;
464  }
465  if (ran > brr[2]) {
466  isw = 4;
467  }
468  *nopt = isw;
469  // **
470  // ** EVAPORATION
471  // **
472  rmnve1 = massPionPlus;
473  rmnve2 = massPionMinus;
474  if (isw == 2) {
475  rmnve1 = massPionZero;
476  }
477  if (isw == 2) {
478  rmnve2 = massPionZero;
479  }
480  if (isw == 3) {
481  rmnve2 = massPionZero;
482  }
483  if (isw == 4) {
484  rmnve1 = massGamma;
485  }
486  if (isw == 4) {
487  rmnve2 = massGamma;
488  }
489  ek = massNeutron + massAntiNeutron - rmnve1 - rmnve2;
490  tkin = ExNu(ek);
491  ek -= tkin;
492  if (ek < G4float(1e-4)) {
493  ek = G4float(1e-4);
494  }
495  ek /= G4float(2.);
496  en = ek + (rmnve1 + rmnve2) / G4float(2.);
497  r__1 = en * en - rmnve1 * rmnve2;
498  pcm = r__1 > 0 ? std::sqrt(r__1) : 0;
499  pv[2].SetZero();
500  pv[2].SetMass( rmnve1 );
501  pv[3].SetZero();
502  pv[3].SetMass( rmnve2 );
503  if (isw > 3) {
504  pv[2].SetMass( 0. );
505  pv[3].SetMass( 0. );
506  }
507  pv[2].SetEnergyAndUpdate( std::sqrt(pv[2].GetMass()*pv[2].GetMass()+pcm*pcm) );
508  pv[2].SetTOF( result.GetTOF() );
509  pv[3].SetEnergy( std::sqrt(pv[3].GetMass()*pv[3].GetMass()+pcm*pcm) );
510  pv[3].SetMomentumAndUpdate( -pv[2].GetMomentum().x(), -pv[2].GetMomentum().y(), -pv[2].GetMomentum().z() );
511  pv[3].SetTOF( result.GetTOF() );
512  switch ((int)isw) {
513  case 1:
516  break;
517  case 2:
520  break;
521  case 3:
524  break;
525  case 4:
528  break;
529  default:
530  break;
531  }
532  nt = 3;
533  if (targetAtomicMass >= G4float(1.5)) {
534  cfa = (targetAtomicMass - G4float(1.)) / G4float(120.) *
535  G4float(.025) * G4Exp(-G4double(targetAtomicMass - G4float(1.)) /
536  G4float(120.));
537  targ = G4float(1.);
538  tex = evapEnergy1;
539  if (tex >= G4float(.001)) {
540  black = (targ * G4float(1.25) +
542  Poisso(black, &nbl);
543  if (G4float(G4int(targ) + nbl) > targetAtomicMass) {
544  nbl = G4int(targetAtomicMass - targ);
545  }
546  if (nt + nbl > (MAX_SECONDARIES - 2)) {
547  nbl = (MAX_SECONDARIES - 2) - nt;
548  }
549  if (nbl > 0) {
550  ekin = tex / nbl;
551  ekin2 = G4float(0.);
552  for (i = 1; i <= nbl; ++i) {
553  if (nt == (MAX_SECONDARIES - 2)) {
554  continue;
555  }
556  if (ekin2 > tex) {
557  break;
558  }
559  ran1 = G4UniformRand();
560  Normal(&ran2);
561  ekin1 = -G4double(ekin) * G4Log(ran1) -
562  cfa * (ran2 * G4float(.5) + G4float(1.));
563  if (ekin1 < G4float(0.)) {
564  ekin1 = G4Log(ran1) * G4float(-.01);
565  }
566  ekin1 *= G4float(1.);
567  ekin2 += ekin1;
568  if (ekin2 > tex) {
569  ekin1 = tex - (ekin2 - ekin1);
570  }
571  if (ekin1 < G4float(0.)) {
572  ekin1 = G4float(.001);
573  }
574  ipa1 = pdefNeutron;
575  pnrat = G4float(1.) - targetCharge / targetAtomicMass;
576  if (G4UniformRand() > pnrat) {
577  ipa1 = pdefProton;
578  }
579  ++nt;
580  pv[nt].SetZero();
581  pv[nt].SetMass( ipa1->GetPDGMass()/GeV );
582  pv[nt].SetKineticEnergyAndUpdate( ekin1 );
583  pv[nt].SetTOF( result.GetTOF() );
584  pv[nt].SetParticleDef( ipa1 );
585  }
586  if (targetAtomicMass >= G4float(230.) && ek <= G4float(2.)) {
587  ii = nt + 1;
588  kk = 0;
589  eka = ek;
590  if (eka > G4float(1.)) {
591  eka *= eka;
592  }
593  if (eka < G4float(.1)) {
594  eka = G4float(.1);
595  }
596  ika = G4int(G4float(3.6) / eka);
597  for (i = 1; i <= nt; ++i) {
598  --ii;
599  if (pv[ii].GetParticleDef() != pdefProton) {
600  continue;
601  }
602  ipa1 = pdefNeutron;
603  pv[ii].SetMass( ipa1->GetPDGMass()/GeV );
604  pv[ii].SetParticleDef( ipa1 );
605  ++kk;
606  if (kk > ika) {
607  break;
608  }
609  }
610  }
611  }
612  }
613  // **
614  // ** THEN ALSO DEUTERONS, TRITONS AND ALPHAS
615  // **
616  tex = evapEnergy3;
617  if (tex >= G4float(.001)) {
618  black = (targ * G4float(1.25) + G4float(1.5)) * evapEnergy3 /
620  Poisso(black, &nbl);
621  if (nt + nbl > (MAX_SECONDARIES - 2)) {
622  nbl = (MAX_SECONDARIES - 2) - nt;
623  }
624  if (nbl > 0) {
625  ekin = tex / nbl;
626  ekin2 = G4float(0.);
627  for (i = 1; i <= nbl; ++i) {
628  if (nt == (MAX_SECONDARIES - 2)) {
629  continue;
630  }
631  if (ekin2 > tex) {
632  break;
633  }
634  ran1 = G4UniformRand();
635  Normal(&ran2);
636  ekin1 = -G4double(ekin) * G4Log(ran1) -
637  cfa * (ran2 * G4float(.5) + G4float(1.));
638  if (ekin1 < G4float(0.)) {
639  ekin1 = G4Log(ran1) * G4float(-.01);
640  }
641  ekin1 *= G4float(1.);
642  ekin2 += ekin1;
643  if (ekin2 > tex) {
644  ekin1 = tex - (ekin2 - ekin1);
645  }
646  if (ekin1 < G4float(0.)) {
647  ekin1 = G4float(.001);
648  }
649  ran = G4UniformRand();
650  inve = pdefDeuteron;
651  if (ran > G4float(.6)) {
652  inve = pdefTriton;
653  }
654  if (ran > G4float(.9)) {
655  inve = pdefAlpha;
656  }
657  ++nt;
658  pv[nt].SetZero();
659  pv[nt].SetMass( inve->GetPDGMass()/GeV );
660  pv[nt].SetKineticEnergyAndUpdate( ekin1 );
661  pv[nt].SetTOF( result.GetTOF() );
662  pv[nt].SetParticleDef( inve );
663  }
664  }
665  }
666  }
667  result = pv[2];
668  if (nt == 2) {
669  return;
670  }
671  for (i = 3; i <= nt; ++i) {
672  if (ntot >= MAX_SECONDARIES) {
673  return;
674  }
675  eve[ntot++] = pv[i];
676  }
677 
678 } // AntiNeutronAnnihilation
void SetParticleDef(G4ParticleDefinition *c)
G4ParticleDefinition * GetParticleDef()
float G4float
Definition: G4Types.hh:77
#define G4ThreadLocal
Definition: tls.hh:89
int G4int
Definition: G4Types.hh:78
void SetKineticEnergyAndUpdate(G4double ekin)
Double_t y
#define G4UniformRand()
Definition: Randomize.hh:97
#define MAX_SECONDARIES
static const double GeV
Definition: G4SIunits.hh:214
TString targ[ntarg]
void SetEnergyAndUpdate(G4double e)
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
void SetMomentumAndUpdate(G4ParticleMomentum mom)
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ AtRestDoIt()

G4VParticleChange * G4AntiNeutronAnnihilationAtRest::AtRestDoIt ( const G4Track &  track,
const G4Step &   
)
virtual

Reimplemented from G4VRestProcess.

Definition at line 151 of file G4AntiNeutronAnnihilationAtRest.cc.

160 {
161 
162 // Initialize ParticleChange
163 // all members of G4VParticleChange are set to equal to
164 // corresponding member in G4Track
165 
166  aParticleChange.Initialize(track);
167 
168 // Store some global quantities that depend on current material and particle
169 
170  globalTime = track.GetGlobalTime()/s;
171  G4Material * aMaterial = track.GetMaterial();
172  const G4int numberOfElements = aMaterial->GetNumberOfElements();
173  const G4ElementVector* theElementVector = aMaterial->GetElementVector();
174 
175  const G4double* theAtomicNumberDensity = aMaterial->GetAtomicNumDensityVector();
176  G4double normalization = 0;
177  for ( G4int i1=0; i1 < numberOfElements; i1++ )
178  {
179  normalization += theAtomicNumberDensity[i1] ; // change when nucleon specific
180  // probabilities are included.
181  }
182  G4double runningSum= 0.;
183  G4double random = G4UniformRand()*normalization;
184  for ( G4int i2=0; i2 < numberOfElements; i2++ )
185  {
186  runningSum += theAtomicNumberDensity[i2]; // change when nucleon specific
187  // probabilities are included.
188  if (random<=runningSum)
189  {
190  targetCharge = G4double( ((*theElementVector)[i2])->GetZ());
191  targetAtomicMass = (*theElementVector)[i2]->GetN();
192  }
193  }
194  if (random>runningSum)
195  {
196  targetCharge = G4double( ((*theElementVector)[numberOfElements-1])->GetZ());
197  targetAtomicMass = (*theElementVector)[numberOfElements-1]->GetN();
198  }
199 
200  if (verboseLevel>1) {
201  G4cout << "G4AntiNeutronAnnihilationAtRest::AtRestDoIt is invoked " <<G4endl;
202  }
203 
204  G4ParticleMomentum momentum;
205  G4float localtime;
206 
207  G4ThreeVector position = track.GetPosition();
208 
209  GenerateSecondaries(); // Generate secondaries
210 
211  aParticleChange.SetNumberOfSecondaries( ngkine );
212 
213  for ( G4int isec = 0; isec < ngkine; isec++ ) {
214  G4DynamicParticle* aNewParticle = new G4DynamicParticle;
215  aNewParticle->SetDefinition( gkin[isec].GetParticleDef() );
216  aNewParticle->SetMomentum( gkin[isec].GetMomentum() * GeV );
217 
218  localtime = globalTime + gkin[isec].GetTOF();
219 
220  G4Track* aNewTrack = new G4Track( aNewParticle, localtime*s, position );
221  aNewTrack->SetTouchableHandle(track.GetTouchableHandle());
222  aParticleChange.AddSecondary( aNewTrack );
223 
224  }
225 
226  aParticleChange.ProposeLocalEnergyDeposit( 0.0*GeV );
227 
228  aParticleChange.ProposeTrackStatus(fStopAndKill); // Kill the incident AntiNeutron
229 
230 // clear InteractionLengthLeft
231 
233 
234  return &aParticleChange;
235 
236 }
void SetMomentum(const G4ThreeVector &momentum)
G4int verboseLevel
Definition: G4VProcess.hh:368
std::vector< G4Element * > G4ElementVector
float G4float
Definition: G4Types.hh:77
virtual void ResetNumberOfInteractionLengthLeft()
Definition: G4VProcess.cc:95
int G4int
Definition: G4Types.hh:78
static const double s
Definition: G4SIunits.hh:168
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:216
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
static const double GeV
Definition: G4SIunits.hh:214
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
#define G4endl
Definition: G4ios.hh:61
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:190
double G4double
Definition: G4Types.hh:76
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ AtRestGetPhysicalInteractionLength()

G4double G4AntiNeutronAnnihilationAtRest::AtRestGetPhysicalInteractionLength ( const G4Track &  track,
G4ForceCondition *  condition 
)
virtual

Reimplemented from G4VRestProcess.

Definition at line 125 of file G4AntiNeutronAnnihilationAtRest.cc.

129 {
130  // beggining of tracking
132 
133  // condition is set to "Not Forced"
134  *condition = NotForced;
135 
136  // get mean life time
138 
139  if ((currentInteractionLength <0.0) || (verboseLevel>2)){
140  G4cout << "G4AntiNeutronAnnihilationAtRestProcess::AtRestGetPhysicalInteractionLength ";
141  G4cout << "[ " << GetProcessName() << "]" <<G4endl;
142  track.GetDynamicParticle()->DumpInfo();
143  G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl;
144  G4cout << "MeanLifeTime = " << currentInteractionLength/ns << "[ns]" <<G4endl;
145  }
146 
148 
149 }
G4double condition(const G4ErrorSymMatrix &m)
G4double GetMeanLifeTime(const G4Track &, G4ForceCondition *)
G4int verboseLevel
Definition: G4VProcess.hh:368
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:293
virtual void ResetNumberOfInteractionLengthLeft()
Definition: G4VProcess.cc:95
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4GLOB_DLL std::ostream G4cout
G4double currentInteractionLength
Definition: G4VProcess.hh:297
#define G4endl
Definition: G4ios.hh:61
#define ns
Definition: xmlparse.cc:614
Here is the call graph for this function:

◆ BuildPhysicsTable()

void G4AntiNeutronAnnihilationAtRest::BuildPhysicsTable ( const G4ParticleDefinition p)
virtual

Reimplemented from G4VProcess.

Definition at line 96 of file G4AntiNeutronAnnihilationAtRest.cc.

97 {
99 }
static G4HadronicProcessStore * Instance()
void PrintInfo(const G4ParticleDefinition *)
Here is the call graph for this function:

◆ ExNu()

G4double G4AntiNeutronAnnihilationAtRest::ExNu ( G4float  ek1)
private

Definition at line 681 of file G4AntiNeutronAnnihilationAtRest.cc.

682 {
683  G4float ret_val, r__1;
684 
685  static G4ThreadLocal G4float cfa, gfa, ran1, ran2, ekin1, atno3;
686  static G4ThreadLocal G4int magic;
687  static G4ThreadLocal G4float fpdiv;
688 
689  // *** NUCLEAR EVAPORATION AS FUNCTION OF ATOMIC NUMBER ATNO ***
690  // *** AND KINETIC ENERGY EKIN OF PRIMARY PARTICLE ***
691  // *** NVE 04-MAR-1988 CERN GENEVA ***
692  // ORIGIN : H.FESEFELDT (10-DEC-1986)
693 
694  ret_val = G4float(0.);
695  if (targetAtomicMass >= G4float(1.5)) {
696  magic = 0;
697  if (G4int(targetCharge + G4float(.1)) == 82) {
698  magic = 1;
699  }
700  ekin1 = ek1;
701  if (ekin1 < G4float(.1)) {
702  ekin1 = G4float(.1);
703  }
704  if (ekin1 > G4float(4.)) {
705  ekin1 = G4float(4.);
706  }
707  // ** 0.35 VALUE AT 1 GEV
708  // ** 0.05 VALUE AT 0.1 GEV
709  cfa = G4float(.13043478260869565);
710  cfa = cfa * G4Log(ekin1) + G4float(.35);
711  if (cfa < G4float(.15)) {
712  cfa = G4float(.15);
713  }
714  ret_val = cfa * G4float(7.716) * G4Exp(-G4double(cfa));
715  atno3 = targetAtomicMass;
716  if (atno3 > G4float(120.)) {
717  atno3 = G4float(120.);
718  }
719  cfa = (atno3 - G4float(1.)) /
720  G4float(120.) * G4Exp(-G4double(atno3 - G4float(1.)) / G4float(120.));
721  ret_val *= cfa;
722  r__1 = ekin1;
723  fpdiv = G4float(1.) - r__1 * r__1 * G4float(.25);
724  if (fpdiv < G4float(.5)) {
725  fpdiv = G4float(.5);
726  }
727  gfa = (targetAtomicMass - G4float(1.)) /
728  G4float(70.) * G4float(2.) *
730  evapEnergy1 = ret_val * fpdiv;
731  evapEnergy3 = ret_val - evapEnergy1;
732  Normal(&ran1);
733  Normal(&ran2);
734  if (magic == 1) {
735  ran1 = G4float(0.);
736  ran2 = G4float(0.);
737  }
738  evapEnergy1 *= ran1 * gfa + G4float(1.);
739  if (evapEnergy1 < G4float(0.)) {
740  evapEnergy1 = G4float(0.);
741  }
742  evapEnergy3 *= ran2 * gfa + G4float(1.);
743  if (evapEnergy3 < G4float(0.)) {
744  evapEnergy3 = G4float(0.);
745  }
746 
747  // Loop checking, 06-Aug-2015, Vladimir Ivanchenko
748  while ((ret_val = evapEnergy1 + evapEnergy3) >= ek1) {
749  evapEnergy1 *= G4float(1.) - G4UniformRand() * G4float(.5);
750  evapEnergy3 *= G4float(1.) - G4UniformRand() * G4float(.5);
751  }
752  }
753  return ret_val;
754 
755 } // ExNu
float G4float
Definition: G4Types.hh:77
#define G4ThreadLocal
Definition: tls.hh:89
int G4int
Definition: G4Types.hh:78
#define G4UniformRand()
Definition: Randomize.hh:97
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GenerateSecondaries()

void G4AntiNeutronAnnihilationAtRest::GenerateSecondaries ( )
private

Definition at line 239 of file G4AntiNeutronAnnihilationAtRest.cc.

240 {
241  static G4ThreadLocal G4int index;
242  static G4ThreadLocal G4int l;
243  static G4ThreadLocal G4int nopt;
244  static G4ThreadLocal G4int i;
245  // DHW 15 May 2011: unused: static G4ParticleDefinition* jnd;
246 
247  for (i = 1; i <= MAX_SECONDARIES; ++i) {
248  pv[i].SetZero();
249  }
250 
251 
252  ngkine = 0; // number of generated secondary particles
253  ntot = 0;
254  result.SetZero();
257  result.SetTOF( 0. );
259 
260  // *** SELECT PROCESS FOR CURRENT PARTICLE ***
261 
263 
264  // *** CHECK WHETHER THERE ARE NEW PARTICLES GENERATED ***
265  if (ntot != 0 || result.GetParticleDef() != pdefAntiNeutron) {
266  // *** CURRENT PARTICLE IS NOT THE SAME AS IN THE BEGINNING OR/AND ***
267  // *** ONE OR MORE SECONDARIES HAVE BEEN GENERATED ***
268 
269  // --- INITIAL PARTICLE TYPE HAS BEEN CHANGED ==> PUT NEW TYPE ON ---
270  // --- THE GEANT TEMPORARY STACK ---
271 
272  // --- PUT PARTICLE ON THE STACK ---
273  gkin[0] = result;
274  gkin[0].SetTOF( result.GetTOF() * 5e-11 );
275  ngkine = 1;
276 
277  // --- ALL QUANTITIES ARE TAKEN FROM THE GHEISHA STACK WHERE THE ---
278  // --- CONVENTION IS THE FOLLOWING ---
279 
280  // --- ONE OR MORE SECONDARIES HAVE BEEN GENERATED ---
281  for (l = 1; l <= ntot; ++l) {
282  index = l - 1;
283  // DHW 15 May 2011: unused: jnd = eve[index].GetParticleDef();
284 
285  // --- ADD PARTICLE TO THE STACK IF STACK NOT YET FULL ---
286  if (ngkine < MAX_SECONDARIES) {
287  gkin[ngkine] = eve[index];
288  gkin[ngkine].SetTOF( eve[index].GetTOF() * 5e-11 );
289  ++ngkine;
290  }
291  }
292  }
293  else {
294  // --- NO SECONDARIES GENERATED AND PARTICLE IS STILL THE SAME ---
295  // --- ==> COPY EVERYTHING BACK IN THE CURRENT GEANT STACK ---
296  ngkine = 0;
297  ntot = 0;
298  globalTime += result.GetTOF() * G4float(5e-11);
299  }
300 
301  // --- LIMIT THE VALUE OF NGKINE IN CASE OF OVERFLOW ---
303 
304 } // GenerateSecondaries
void SetParticleDef(G4ParticleDefinition *c)
Int_t index
G4ParticleDefinition * GetParticleDef()
float G4float
Definition: G4Types.hh:77
#define G4ThreadLocal
Definition: tls.hh:89
int G4int
Definition: G4Types.hh:78
void SetKineticEnergyAndUpdate(G4double ekin)
#define MAX_SECONDARIES
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetMeanLifeTime()

G4double G4AntiNeutronAnnihilationAtRest::GetMeanLifeTime ( const G4Track &  ,
G4ForceCondition *   
)
inlinevirtual

Implements G4VRestProcess.

Definition at line 72 of file G4AntiNeutronAnnihilationAtRest.hh.

73  {return 0.0;}
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetNumberOfSecondaries()

G4int G4AntiNeutronAnnihilationAtRest::GetNumberOfSecondaries ( )

Definition at line 112 of file G4AntiNeutronAnnihilationAtRest.cc.

113 {
114  return ( ngkine );
115 
116 }
Here is the caller graph for this function:

◆ GetSecondaryKinematics()

G4GHEKinematicsVector * G4AntiNeutronAnnihilationAtRest::GetSecondaryKinematics ( )

Definition at line 119 of file G4AntiNeutronAnnihilationAtRest.cc.

120 {
121  return ( &gkin[0] );
122 
123 }
Here is the caller graph for this function:

◆ IsApplicable()

G4bool G4AntiNeutronAnnihilationAtRest::IsApplicable ( const G4ParticleDefinition particle)
virtual

Reimplemented from G4VProcess.

Definition at line 103 of file G4AntiNeutronAnnihilationAtRest.cc.

106 {
107  return ( &particle == pdefAntiNeutron );
108 
109 }

◆ NFac()

G4int G4AntiNeutronAnnihilationAtRest::NFac ( G4int  n)
private

Definition at line 382 of file G4AntiNeutronAnnihilationAtRest.cc.

383 {
384  G4int ret_val;
385 
386  static G4ThreadLocal G4int i, j;
387 
388  // *** NVE 16-MAR-1988 CERN GENEVA ***
389  // ORIGIN : H.FESEFELDT (27-OCT-1983)
390 
391  ret_val = 1;
392  j = n;
393  if (j > 1) {
394  if (j > 10) {
395  j = 10;
396  }
397  for (i = 2; i <= j; ++i) {
398  ret_val *= i;
399  }
400  }
401  return ret_val;
402 
403 } // NFac
#define G4ThreadLocal
Definition: tls.hh:89
int G4int
Definition: G4Types.hh:78
Char_t n[5]
Here is the caller graph for this function:

◆ Normal()

void G4AntiNeutronAnnihilationAtRest::Normal ( G4float ran)
private

Definition at line 406 of file G4AntiNeutronAnnihilationAtRest.cc.

407 {
408  static G4ThreadLocal G4int i;
409 
410  // *** NVE 14-APR-1988 CERN GENEVA ***
411  // ORIGIN : H.FESEFELDT (27-OCT-1983)
412 
413  *ran = G4float(-6.);
414  for (i = 1; i <= 12; ++i) {
415  *ran += G4UniformRand();
416  }
417 
418 } // Normal
float G4float
Definition: G4Types.hh:77
#define G4ThreadLocal
Definition: tls.hh:89
int G4int
Definition: G4Types.hh:78
#define G4UniformRand()
Definition: Randomize.hh:97
Here is the caller graph for this function:

◆ operator=()

G4AntiNeutronAnnihilationAtRest& G4AntiNeutronAnnihilationAtRest::operator= ( const G4AntiNeutronAnnihilationAtRest right)
private

◆ Poisso()

void G4AntiNeutronAnnihilationAtRest::Poisso ( G4float  xav,
G4int iran 
)
private

Definition at line 307 of file G4AntiNeutronAnnihilationAtRest.cc.

308 {
309  static G4ThreadLocal G4int i;
310  static G4ThreadLocal G4float r, p1, p2, p3;
311  static G4ThreadLocal G4int fivex;
312  static G4ThreadLocal G4float rr, ran, rrr, ran1;
313 
314  // *** GENERATION OF POISSON DISTRIBUTION ***
315  // *** NVE 16-MAR-1988 CERN GENEVA ***
316  // ORIGIN : H.FESEFELDT (27-OCT-1983)
317 
318  // --- USE NORMAL DISTRIBUTION FOR <X> > 9.9 ---
319  if (xav > G4float(9.9)) {
320  // ** NORMAL DISTRIBUTION WITH SIGMA**2 = <X>
321  Normal(&ran1);
322  ran1 = xav + ran1 * std::sqrt(xav);
323  *iran = G4int(ran1);
324  if (*iran < 0) {
325  *iran = 0;
326  }
327  }
328  else {
329  fivex = G4int(xav * G4float(5.));
330  *iran = 0;
331  if (fivex > 0) {
332  r = G4Exp(-G4double(xav));
333  ran1 = G4UniformRand();
334  if (ran1 > r) {
335  rr = r;
336  for (i = 1; i <= fivex; ++i) {
337  ++(*iran);
338  if (i <= 5) {
339  rrr = G4Pow::GetInstance()->powN(xav, i) / NFac(i);
340  }
341  // ** STIRLING' S FORMULA FOR LARGE NUMBERS
342  if (i > 5) {
343  rrr = G4Exp(i * G4Log(xav) -
344  (i + G4float(.5)) * G4Log(i * G4float(1.)) +
345  i - G4float(.9189385));
346  }
347  rr += r * rrr;
348  if (ran1 <= rr) {
349  break;
350  }
351  }
352  }
353  }
354  else {
355  // ** FOR VERY SMALL XAV TRY IRAN=1,2,3
356  p1 = xav * G4Exp(-G4double(xav));
357  p2 = xav * p1 / G4float(2.);
358  p3 = xav * p2 / G4float(3.);
359  ran = G4UniformRand();
360  if (ran >= p3) {
361  if (ran >= p2) {
362  if (ran >= p1) {
363  *iran = 0;
364  }
365  else {
366  *iran = 1;
367  }
368  }
369  else {
370  *iran = 2;
371  }
372  }
373  else {
374  *iran = 3;
375  }
376  }
377  }
378 
379 } // Poisso
static G4Pow * GetInstance()
Definition: G4Pow.cc:55
G4double powN(G4double x, G4int n) const
Definition: G4Pow.cc:128
float G4float
Definition: G4Types.hh:77
#define G4ThreadLocal
Definition: tls.hh:89
int G4int
Definition: G4Types.hh:78
#define G4UniformRand()
Definition: Randomize.hh:97
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ PreparePhysicsTable()

void G4AntiNeutronAnnihilationAtRest::PreparePhysicsTable ( const G4ParticleDefinition p)
virtual

Reimplemented from G4VProcess.

Definition at line 91 of file G4AntiNeutronAnnihilationAtRest.cc.

92 {
94 }
static G4HadronicProcessStore * Instance()
void RegisterParticleForExtraProcess(G4VProcess *, const G4ParticleDefinition *)
Here is the call graph for this function:

Member Data Documentation

◆ evapEnergy1

G4float G4AntiNeutronAnnihilationAtRest::evapEnergy1
private

Definition at line 107 of file G4AntiNeutronAnnihilationAtRest.hh.

◆ evapEnergy3

G4float G4AntiNeutronAnnihilationAtRest::evapEnergy3
private

Definition at line 108 of file G4AntiNeutronAnnihilationAtRest.hh.

◆ eve

G4GHEKinematicsVector* G4AntiNeutronAnnihilationAtRest::eve
private

Definition at line 104 of file G4AntiNeutronAnnihilationAtRest.hh.

◆ gkin

G4GHEKinematicsVector* G4AntiNeutronAnnihilationAtRest::gkin
private

Definition at line 105 of file G4AntiNeutronAnnihilationAtRest.hh.

◆ globalTime

G4float G4AntiNeutronAnnihilationAtRest::globalTime
private

Definition at line 95 of file G4AntiNeutronAnnihilationAtRest.hh.

◆ massAntiNeutron

G4float G4AntiNeutronAnnihilationAtRest::massAntiNeutron
private

Definition at line 119 of file G4AntiNeutronAnnihilationAtRest.hh.

◆ massGamma

G4float G4AntiNeutronAnnihilationAtRest::massGamma
private

Definition at line 118 of file G4AntiNeutronAnnihilationAtRest.hh.

◆ massNeutron

G4float G4AntiNeutronAnnihilationAtRest::massNeutron
private

Definition at line 120 of file G4AntiNeutronAnnihilationAtRest.hh.

◆ massPionMinus

G4float G4AntiNeutronAnnihilationAtRest::massPionMinus
private

Definition at line 115 of file G4AntiNeutronAnnihilationAtRest.hh.

◆ massPionPlus

G4float G4AntiNeutronAnnihilationAtRest::massPionPlus
private

Definition at line 117 of file G4AntiNeutronAnnihilationAtRest.hh.

◆ massPionZero

G4float G4AntiNeutronAnnihilationAtRest::massPionZero
private

Definition at line 116 of file G4AntiNeutronAnnihilationAtRest.hh.

◆ ngkine

G4int G4AntiNeutronAnnihilationAtRest::ngkine
private

Definition at line 110 of file G4AntiNeutronAnnihilationAtRest.hh.

◆ ntot

G4int G4AntiNeutronAnnihilationAtRest::ntot
private

Definition at line 112 of file G4AntiNeutronAnnihilationAtRest.hh.

◆ pdefAlpha

G4ParticleDefinition* G4AntiNeutronAnnihilationAtRest::pdefAlpha
private

Definition at line 131 of file G4AntiNeutronAnnihilationAtRest.hh.

◆ pdefAntiNeutron

G4ParticleDefinition* G4AntiNeutronAnnihilationAtRest::pdefAntiNeutron
private

Definition at line 128 of file G4AntiNeutronAnnihilationAtRest.hh.

◆ pdefDeuteron

G4ParticleDefinition* G4AntiNeutronAnnihilationAtRest::pdefDeuteron
private

Definition at line 129 of file G4AntiNeutronAnnihilationAtRest.hh.

◆ pdefGamma

G4ParticleDefinition* G4AntiNeutronAnnihilationAtRest::pdefGamma
private

Definition at line 122 of file G4AntiNeutronAnnihilationAtRest.hh.

◆ pdefNeutron

G4ParticleDefinition* G4AntiNeutronAnnihilationAtRest::pdefNeutron
private

Definition at line 127 of file G4AntiNeutronAnnihilationAtRest.hh.

◆ pdefPionMinus

G4ParticleDefinition* G4AntiNeutronAnnihilationAtRest::pdefPionMinus
private

Definition at line 125 of file G4AntiNeutronAnnihilationAtRest.hh.

◆ pdefPionPlus

G4ParticleDefinition* G4AntiNeutronAnnihilationAtRest::pdefPionPlus
private

Definition at line 123 of file G4AntiNeutronAnnihilationAtRest.hh.

◆ pdefPionZero

G4ParticleDefinition* G4AntiNeutronAnnihilationAtRest::pdefPionZero
private

Definition at line 124 of file G4AntiNeutronAnnihilationAtRest.hh.

◆ pdefProton

G4ParticleDefinition* G4AntiNeutronAnnihilationAtRest::pdefProton
private

Definition at line 126 of file G4AntiNeutronAnnihilationAtRest.hh.

◆ pdefTriton

G4ParticleDefinition* G4AntiNeutronAnnihilationAtRest::pdefTriton
private

Definition at line 130 of file G4AntiNeutronAnnihilationAtRest.hh.

◆ pv

G4GHEKinematicsVector* G4AntiNeutronAnnihilationAtRest::pv
private

Definition at line 103 of file G4AntiNeutronAnnihilationAtRest.hh.

◆ result

G4GHEKinematicsVector G4AntiNeutronAnnihilationAtRest::result
private

Definition at line 113 of file G4AntiNeutronAnnihilationAtRest.hh.

◆ targetAtomicMass

G4float G4AntiNeutronAnnihilationAtRest::targetAtomicMass
private

Definition at line 98 of file G4AntiNeutronAnnihilationAtRest.hh.

◆ targetCharge

G4float G4AntiNeutronAnnihilationAtRest::targetCharge
private

Definition at line 101 of file G4AntiNeutronAnnihilationAtRest.hh.


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