Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4QuasiElasticChannel.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// $Id$
28//
29
30// Author : Gunter Folger March 2007
31// Modified by Mikhail Kossov. Apr2009, E/M conservation: ResidualNucleus is added (ResNuc)
32// Class Description
33// Final state production model for theoretical models of hadron inelastic
34// quasi elastic scattering in geant4;
35// Class Description - End
36//
37// Modified:
38// 20110805 M. Kelsey -- Follow change to G4V3DNucleus::GetNucleons()
39// 20110808 M. Kelsey -- Move #includes from .hh, add many missing
40
42
43#include "G4Fancy3DNucleus.hh"
44#include "G4DynamicParticle.hh"
45#include "G4HadTmpUtil.hh" /* lrint */
46#include "G4KineticTrack.hh"
48#include "G4LorentzVector.hh"
49#include "G4Neutron.hh"
50#include "G4Nucleon.hh"
51#include "G4Nucleus.hh"
53#include "G4ParticleTable.hh"
54#include "G4QuasiElRatios.hh"
55#include "globals.hh"
56#include <vector>
57
58//#define debug_scatter
59
60
62 : theQuasiElastic(G4QuasiElRatios::GetPointer()),
63 the3DNucleus(new G4Fancy3DNucleus) {}
64
66{
67 delete the3DNucleus;
68}
69
71 const G4DynamicParticle & thePrimary)
72{
73 #ifdef debug_scatter
74 G4cout << "G4QuasiElasticChannel:: P=" << thePrimary.GetTotalMomentum()
75 << ", pPDG=" << thePrimary.GetDefinition()->GetPDGEncoding()
76 << ", Z = " << theNucleus.GetZ_asInt())
77 << ", N = " << theNucleus.GetN_asInt())
78 << ", A = " << theNucleus.GetA_asInt() << G4endl;
79 #endif
80
81 std::pair<G4double,G4double> ratios;
82 ratios=theQuasiElastic->GetRatios(thePrimary.GetTotalMomentum(),
83 thePrimary.GetDefinition()->GetPDGEncoding(),
84 theNucleus.GetZ_asInt(),
85 theNucleus.GetN_asInt());
86 #ifdef debug_scatter
87 G4cout << "G4QuasiElasticChannel::ratios " << ratios.first << " x " <<ratios.second
88 << " = " << ratios.first*ratios.second << G4endl;
89 #endif
90
91 return ratios.first*ratios.second;
92}
93
95 const G4DynamicParticle & thePrimary)
96{
97 G4int A=theNucleus.GetA_asInt();
98 G4int Z=theNucleus.GetZ_asInt();
99 // build Nucleus and choose random nucleon to scatter with
100 the3DNucleus->Init(theNucleus.GetA_asInt(),theNucleus.GetZ_asInt());
101 const std::vector<G4Nucleon>& nucleons=the3DNucleus->GetNucleons();
102 G4double targetNucleusMass=the3DNucleus->GetMass();
103 G4LorentzVector targetNucleus4Mom(0.,0.,0.,targetNucleusMass);
104 G4int index;
105 do {
106 index=G4lrint((A-1)*G4UniformRand());
107 } while (index < 0 || index >= static_cast<G4int>(nucleons.size()));
108
109 G4ParticleDefinition * pDef= nucleons[index].GetDefinition();
110
111 G4int resA=A - 1;
112 G4int resZ=Z - static_cast<int>(pDef->GetPDGCharge());
113 G4ParticleDefinition* resDef;
114 G4double residualNucleusMass;
115 if(resZ)
116 {
117 resDef=G4ParticleTable::GetParticleTable()->FindIon(resZ,resA,0,resZ);
118 residualNucleusMass=resDef->GetPDGMass();
119 }
120 else {
121 resDef=G4Neutron::Neutron();
122 residualNucleusMass=resA * G4Neutron::Neutron()->GetPDGMass();
123 }
124 #ifdef debug_scatter
125 G4cout<<"G4QElChan::Scatter: neutron - proton? A ="<<A<<", Z="<<Z<<", projName="
126 <<pDef->GetParticleName()<<G4endl;
127 #endif
128
129 G4LorentzVector pNucleon=nucleons[index].Get4Momentum();
130 G4double residualNucleusEnergy=std::sqrt(sqr(residualNucleusMass) +
131 pNucleon.vect().mag2());
132 pNucleon.setE(targetNucleusMass-residualNucleusEnergy);
133 G4LorentzVector residualNucleus4Mom=targetNucleus4Mom-pNucleon;
134
135 std::pair<G4LorentzVector,G4LorentzVector> result;
136
137 result=theQuasiElastic->Scatter(pDef->GetPDGEncoding(),pNucleon,
138 thePrimary.GetDefinition()->GetPDGEncoding(),
139 thePrimary.Get4Momentum());
140 G4LorentzVector scatteredHadron4Mom;
141 if (result.first.e() > 0.)
142 scatteredHadron4Mom=result.second;
143 else { //scatter failed
144 //G4cout << "Warning - G4QuasiElasticChannel::Scatter no scattering" << G4endl;
145 //return 0; //no scatter
146 scatteredHadron4Mom=thePrimary.Get4Momentum();
147 residualNucleus4Mom=G4LorentzVector(0.,0.,0.,targetNucleusMass);
149 }
150
151#ifdef debug_scatter
152 G4LorentzVector EpConservation=pNucleon+thePrimary.Get4Momentum()
153 - result.first - result.second;
154 if ( (EpConservation.vect().mag2() > .01*MeV*MeV )
155 || (std::abs(EpConservation.e()) > 0.1 * MeV ) )
156 {
157 G4cout << "Warning - G4QuasiElasticChannel::Scatter E-p non conservation : "
158 << EpConservation << G4endl;
159 }
160#endif
161
163 G4KineticTrack * sPrim=new G4KineticTrack(thePrimary.GetDefinition(),
164 0.,G4ThreeVector(0), scatteredHadron4Mom);
165 ktv->push_back(sPrim);
166 if (result.first.e() > 0.)
167 {
168 G4KineticTrack * sNuc=new G4KineticTrack(pDef, 0.,G4ThreeVector(0), result.first);
169 ktv->push_back(sNuc);
170 }
171 if(resZ || resA==1) // For the only neutron or for tnuclei with Z>0
172 {
173 G4KineticTrack * rNuc=new G4KineticTrack(resDef,
174 0.,G4ThreeVector(0), residualNucleus4Mom);
175 ktv->push_back(rNuc);
176 }
177 else // The residual nucleus consists of only neutrons
178 {
179 residualNucleus4Mom/=resA; // Split 4-mom of A*n system equally
180 for(G4int in=0; in<resA; in++) // Loop over neutrons in A*n system.
181 {
182 G4KineticTrack* rNuc=new G4KineticTrack(resDef,
183 0.,G4ThreeVector(0), residualNucleus4Mom);
184 ktv->push_back(rNuc);
185 }
186 }
187#ifdef debug_scatter
188 G4cout<<"G4QElC::Scat: Nucleon: "<<result.first <<" mass "<<result.first.mag() << G4endl;
189 G4cout<<"G4QElC::Scat: Project: "<<result.second<<" mass "<<result.second.mag()<< G4endl;
190#endif
191 return ktv;
192}
CLHEP::HepLorentzVector G4LorentzVector
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:53
double mag2() const
Hep3Vector vect() const
G4ParticleDefinition * GetDefinition() const
G4LorentzVector Get4Momentum() const
G4double GetTotalMomentum() const
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
G4int GetN_asInt() const
Definition: G4Nucleus.hh:112
G4double GetPDGCharge() const
const G4String & GetParticleName() const
G4ParticleDefinition * FindIon(G4int atomicNumber, G4int atomicMass, G4double excitationEnergy)
static G4ParticleTable * GetParticleTable()
std::pair< G4double, G4double > GetRatios(G4double pIU, G4int prPDG, G4int tgZ, G4int tgN)
std::pair< G4LorentzVector, G4LorentzVector > Scatter(G4int NPDG, G4LorentzVector N4M, G4int pPDG, G4LorentzVector p4M)
G4KineticTrackVector * Scatter(G4Nucleus &theNucleus, const G4DynamicParticle &thePrimary)
G4double GetFraction(G4Nucleus &theNucleus, const G4DynamicParticle &thePrimary)
virtual void Init(G4int theA, G4int theZ)=0
virtual G4double GetMass()=0
virtual const std::vector< G4Nucleon > & GetNucleons()=0
int G4lrint(double ad)
Definition: templates.hh:163
T sqr(const T &x)
Definition: templates.hh:145