Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PiMinusStopMaterial Class Referenceabstract

#include <G4PiMinusStopMaterial.hh>

+ Inheritance diagram for G4PiMinusStopMaterial:

Public Member Functions

 G4PiMinusStopMaterial ()
 
virtual ~G4PiMinusStopMaterial ()
 
virtual std::vector< G4ParticleDefinition * > * DefinitionVector ()
 
virtual std::vector< G4LorentzVector * > * P4Vector (const G4double binding, const G4double mass)
 
virtual G4double FinalNucleons ()=0
 

Protected Member Functions

G4double GenerateAngle (G4double range)
 
G4LorentzVector MakeP4 (G4double p, G4double theta, G4double phi, G4double e)
 
G4double RecoilEnergy (const G4double mass)
 

Protected Attributes

std::vector< G4ParticleDefinition * > * _definitions
 
std::vector< G4LorentzVector * > * _momenta
 
G4DistributionGenerator_distributionE
 
G4DistributionGenerator_distributionAngle
 
G4double theR
 

Detailed Description

Definition at line 49 of file G4PiMinusStopMaterial.hh.

Constructor & Destructor Documentation

◆ G4PiMinusStopMaterial()

G4PiMinusStopMaterial::G4PiMinusStopMaterial ( )

Definition at line 55 of file G4PiMinusStopMaterial.cc.

56{
57 _definitions = 0;
58 _momenta = 0;
61 theR = 0.5;
62}
std::vector< G4ParticleDefinition * > * _definitions
G4DistributionGenerator * _distributionE
std::vector< G4LorentzVector * > * _momenta
G4DistributionGenerator * _distributionAngle

◆ ~G4PiMinusStopMaterial()

G4PiMinusStopMaterial::~G4PiMinusStopMaterial ( )
virtual

Definition at line 67 of file G4PiMinusStopMaterial.cc.

68{
69 if (_definitions != 0) delete _definitions;
70 _definitions = 0;
71
72 //A.R. 26-Jul-2012 Coverity fix
73 if (_momenta != 0) {
74 for (unsigned int i=0; i<_momenta->size(); i++) delete(*_momenta)[i];
75 delete _momenta;
76 }
77
78 delete _distributionE;
79 delete _distributionAngle;
80}

Member Function Documentation

◆ DefinitionVector()

std::vector< G4ParticleDefinition * > * G4PiMinusStopMaterial::DefinitionVector ( )
virtual

Definition at line 82 of file G4PiMinusStopMaterial.cc.

83{
84
85 _definitions->push_back(G4Neutron::Neutron());
86
87 G4double ranflat = G4UniformRand();
88 if (ranflat < theR)
89 { _definitions->push_back(G4Proton::Proton()); }
90 else
91 { _definitions->push_back(G4Neutron::Neutron()); }
92
93 return _definitions;
94
95}
double G4double
Definition: G4Types.hh:64
#define G4UniformRand()
Definition: Randomize.hh:53
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
static G4Proton * Proton()
Definition: G4Proton.cc:93

Referenced by G4PiMinusStopAbsorption::DoAbsorption().

◆ FinalNucleons()

◆ GenerateAngle()

G4double G4PiMinusStopMaterial::GenerateAngle ( G4double  range)
protected

Definition at line 165 of file G4PiMinusStopMaterial.cc.

166{
167 G4double ranflat = G4UniformRand();
168 G4double value = ranflat * x;
169 return value;
170}

Referenced by P4Vector().

◆ MakeP4()

G4LorentzVector G4PiMinusStopMaterial::MakeP4 ( G4double  p,
G4double  theta,
G4double  phi,
G4double  e 
)
protected

Definition at line 172 of file G4PiMinusStopMaterial.cc.

173{
174 // G4LorentzVector p4;
175 G4double px = p * std::sin(theta) * std::cos(phi);
176 G4double py = p * std::sin(theta) * std::sin(phi);
177 G4double pz = p * std::cos(theta);
178 G4LorentzVector p4(px,py,pz,e);
179 return p4;
180}

Referenced by P4Vector().

◆ P4Vector()

std::vector< G4LorentzVector * > * G4PiMinusStopMaterial::P4Vector ( const G4double  binding,
const G4double  mass 
)
virtual

Definition at line 98 of file G4PiMinusStopMaterial.cc.

100{
101 // Generate energy of direct absorption products according to experimental
102 // data. The energy distribution of the two nucleons is assumed to be the
103 // same for protons and neutrons.
104
105 G4double eKin1;
106 G4double eKin2;
107 G4double eRecoil;
108
109 // Assume absorption on two nucleons
110 G4int nNucleons = 2;
111 G4double availableE = G4PionMinus::PionMinus()->GetPDGMass() - nNucleons * binding;
114
115 do
116 {
117 G4double ranflat;
118 G4double p;
119 G4double energy;
120 G4double mass;
121
122 ranflat = G4UniformRand();
123 eKin1 = _distributionE->Generate(ranflat);
124 mass = (*_definitions)[0]->GetPDGMass();
125 energy = eKin1 + mass;
126 p = std::sqrt(energy*energy - mass*mass);
127 G4double theta1 = pi*G4UniformRand();
128 G4double phi1 = GenerateAngle(2.*pi);
129 p1 = MakeP4(p,theta1,phi1,energy);
130
131 ranflat = G4UniformRand();
132 eKin2 = _distributionE->Generate(ranflat);
133 mass = (*_definitions)[1]->GetPDGMass();
134 energy = eKin2 + mass;
135 p = std::sqrt(energy*energy - mass*mass);
136 ranflat = G4UniformRand();
137 G4double opAngle = _distributionAngle->Generate(ranflat);
138 G4double theta2 = theta1 + opAngle;
139 G4double phi2 = phi1 + opAngle;
140
141 p2 = MakeP4(p,theta2,phi2,energy);
142
143 G4double pNucleus = (p1.vect() + p2.vect()).mag();
144 eRecoil = std::sqrt(pNucleus*pNucleus + massNucleus*massNucleus) - massNucleus;
145
146 // ---- Debug
147 // G4cout << " ---- binding = " << binding << ", nucleus mass = " << massNucleus
148 // << ", p nucleus = " << pNucleus << G4endl;
149 // G4cout << "eKin1,2 " << eKin1 << " " << eKin2 << " eRecoil " << eRecoil
150 // << " availableE " << availableE << G4endl;
151 // ----
152
153 } while ((eKin1 + eKin2 + eRecoil) > availableE);
154
155 //A.R. 26-Jul-2012 Coverity fix
156 if (_momenta != 0) {
157 _momenta->push_back(new G4LorentzVector(p1));
158 _momenta->push_back(new G4LorentzVector(p2));
159 }
160
161 return _momenta;
162
163}
CLHEP::HepLorentzVector G4LorentzVector
int G4int
Definition: G4Types.hh:66
Hep3Vector vect() const
G4double Generate(G4double ranflat)
G4double GenerateAngle(G4double range)
G4LorentzVector MakeP4(G4double p, G4double theta, G4double phi, G4double e)
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:98
const G4double pi

Referenced by G4PiMinusStopAbsorption::DoAbsorption().

◆ RecoilEnergy()

G4double G4PiMinusStopMaterial::RecoilEnergy ( const G4double  mass)
protected

Definition at line 182 of file G4PiMinusStopMaterial.cc.

183{
184 G4ThreeVector p(0.,0.,0.);
185
186 //A.R. 26-Jul-2012 Coverity fix
187 if (_momenta != 0) {
188 for (unsigned int i = 0; i< _momenta->size(); i++)
189 {
190 p = p + (*_momenta)[i]->vect();
191 }
192 }
193
194 G4double pNucleus = p.mag();
195 G4double eNucleus = std::sqrt(pNucleus*pNucleus + mass*mass);
196
197 return eNucleus;
198}

Member Data Documentation

◆ _definitions

◆ _distributionAngle

◆ _distributionE

◆ _momenta

◆ theR


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