Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4TritonDecay.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////////////////////////////////////////////////////////////////////////////////
27// //
28// File: G4TritonDecay.cc based on G4AlphaDecay //
29// Author: A.Mereagglia (CENBG), imported in G4repository by L.Desorgher //
30// Date: 24 June 2019, imported in Geant4 on 31 October 2019 //
31// Description: //
32// //
33// //
34////////////////////////////////////////////////////////////////////////////////
35
36#include "G4TritonDecay.hh"
37#include "G4IonTable.hh"
38#include "Randomize.hh"
39#include "G4ThreeVector.hh"
40#include "G4DynamicParticle.hh"
41#include "G4DecayProducts.hh"
43#include "G4SystemOfUnits.hh"
44#include <iostream>
45#include <iomanip>
46
48 const G4double& branch, const G4double& Qvalue,
49 const G4double& excitationE,
50 const G4Ions::G4FloatLevelBase& flb)
51 : G4NuclearDecay("triton decay", Triton, excitationE, flb), transitionQ(Qvalue)
52{
53 SetParent(theParentNucleus); // Store name of parent nucleus, delete G4MT_parent
54 SetBR(branch);
55
57 G4IonTable* theIonTable =
59 G4int daughterZ = theParentNucleus->GetAtomicNumber() - 1;
60 G4int daughterA = theParentNucleus->GetAtomicMass() - 3;
61 SetDaughter(0, theIonTable->GetIon(daughterZ, daughterA, excitationE, flb) );
62 SetDaughter(1, "triton");
63}
64
65
67{}
68
69
71{
72 // Fill G4MT_parent with theParentNucleus (stored by SetParent in ctor)
74
75 // Fill G4MT_daughters with triton and residual nucleus (stored by SetDaughter)
77
78 G4double tritonMass = G4MT_daughters[1]->GetPDGMass();
79 // Excitation energy included in PDG mass
80 G4double nucleusMass = G4MT_daughters[0]->GetPDGMass();
81
82 // Q value was calculated from atomic masses.
83 // Use it to get correct triton energy.
84 G4double cmMomentum = std::sqrt(transitionQ*(transitionQ + 2.*tritonMass)*
85 (transitionQ + 2.*nucleusMass)*
86 (transitionQ + 2.*tritonMass + 2.*nucleusMass) )/
87 (transitionQ + tritonMass + nucleusMass)/2.;
88
89 // Set up final state
90 // parentParticle is set at rest here because boost with correct momentum
91 // is done later
92 G4DynamicParticle parentParticle(G4MT_parent, G4ThreeVector(0,0,0), 0.0);
93 G4DecayProducts* products = new G4DecayProducts(parentParticle);
94
95 G4double costheta = 2.*G4UniformRand()-1.0;
96 G4double sintheta = std::sqrt(1.0 - costheta*costheta);
97 G4double phi = twopi*G4UniformRand()*rad;
98 G4ThreeVector direction(sintheta*std::cos(phi),sintheta*std::sin(phi),
99 costheta);
100
101 G4double KE = std::sqrt(cmMomentum*cmMomentum + tritonMass*tritonMass)
102 - tritonMass;
103 G4DynamicParticle* daughterparticle =
104 new G4DynamicParticle(G4MT_daughters[1], direction, KE, tritonMass);
105 products->PushProducts(daughterparticle);
106
107 KE = std::sqrt(cmMomentum*cmMomentum + nucleusMass*nucleusMass) - nucleusMass;
108 daughterparticle =
109 new G4DynamicParticle(G4MT_daughters[0], -1.0*direction, KE, nucleusMass);
110 products->PushProducts(daughterparticle);
111
112 // Energy conservation check
113 // For triton decays, do final energy check against reaction Q value
114 // which is well-measured using atomic mass differences. Nuclear masses
115 // should not be used since they are not usually directly measured and we
116 // always decay atoms and not fully stripped nuclei.
117 /*
118 G4int nProd = products->entries();
119 G4DynamicParticle* temp = 0;
120 G4double Esum = 0.0;
121 for (G4int i = 0; i < nProd; i++) {
122 temp = products->operator[](i);
123 Esum += temp->GetKineticEnergy();
124 }
125 G4double eCons = (transitionQ - Esum)/keV;
126 if (eCons > 1.e-07) G4cout << " Triton decay check: Ediff (keV) = " << eCons << G4endl;
127 */
128 return products;
129}
130
131
133{
134 G4cout << " G4TritonDecay for parent nucleus " << GetParentName() << G4endl;
135 G4cout << " decays to " << GetDaughterName(0) << " + " << GetDaughterName(1)
136 << " with branching ratio " << GetBR() << "% and Q value "
137 << transitionQ << G4endl;
138}
139
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
G4int PushProducts(G4DynamicParticle *aParticle)
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:522
G4FloatLevelBase
Definition: G4Ions.hh:83
G4int GetAtomicNumber() const
G4int GetAtomicMass() const
G4IonTable * GetIonTable() const
static G4ParticleTable * GetParticleTable()
G4TritonDecay(const G4ParticleDefinition *theParentNucleus, const G4double &theBR, const G4double &Qvalue, const G4double &excitation, const G4Ions::G4FloatLevelBase &flb)
virtual void DumpNuclearInfo()
virtual ~G4TritonDecay()
virtual G4DecayProducts * DecayIt(G4double)
G4ParticleDefinition ** G4MT_daughters
G4double GetBR() const
const G4String & GetParentName() const
void SetBR(G4double value)
void SetNumberOfDaughters(G4int value)
G4ParticleDefinition * G4MT_parent
void CheckAndFillDaughters()
void SetDaughter(G4int anIndex, const G4ParticleDefinition *particle_type)
const G4String & GetDaughterName(G4int anIndex) const
void SetParent(const G4ParticleDefinition *particle_type)