Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4LENDElastic.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 #include "G4LENDElastic.hh"
28 #include "G4PhysicalConstants.hh"
29 #include "G4SystemOfUnits.hh"
30 #include "G4Nucleus.hh"
31 #include "G4ParticleTable.hh"
32 
34 {
35 
36  G4double temp = aTrack.GetMaterial()->GetTemperature();
37 
38  //G4int iZ = int ( aTarg.GetZ() );
39  //G4int iA = int ( aTarg.GetN() );
40  //migrate to integer A and Z (GetN_asInt returns number of neutrons in the nucleus since this)
41  G4int iZ = aTarg.GetZ_asInt();
42  G4int iA = aTarg.GetA_asInt();
43  G4int iM = 0;
44  if ( aTarg.GetIsotope() != NULL ) {
45  iM = aTarg.GetIsotope()->Getm();
46  }
47 
48  G4double ke = aTrack.GetKineticEnergy();
49 
50  //G4HadFinalState* theResult = new G4HadFinalState();
51  G4HadFinalState* theResult = &theParticleChange;
52  theResult->Clear();
53 
54  G4GIDI_target* aTarget = usedTarget_map.find( lend_manager->GetNucleusEncoding( iZ , iA , iM ) )->second->GetTarget();
55  G4double aMu = aTarget->getElasticFinalState( ke*MeV, temp, NULL, NULL );
56 
58  G4double theta = std::acos( aMu );
59  //G4double sinth = std::sin( theta );
60 
61  G4ReactionProduct theNeutron( const_cast<G4ParticleDefinition *>( aTrack.GetDefinition() ) );
62  theNeutron.SetMomentum( aTrack.Get4Momentum().vect() );
63  theNeutron.SetKineticEnergy( ke );
64 
65 //G4cout << "iZ " << iZ << " iA " << iA << G4endl;
66 
67  G4ReactionProduct theTarget( G4ParticleTable::GetParticleTable()->FindIon( iZ , iA , 0 , iZ ) );
68 
69  G4double mass = G4ParticleTable::GetParticleTable()->FindIon( iZ , iA , 0 , iZ )->GetPDGMass();
70 
71 // add Thermal motion
72  G4double kT = k_Boltzmann*temp;
73  G4ThreeVector v ( G4RandGauss::shoot() * std::sqrt( kT*mass )
74  , G4RandGauss::shoot() * std::sqrt( kT*mass )
75  , G4RandGauss::shoot() * std::sqrt( kT*mass ) );
76  theTarget.SetMomentum( v );
77 
78  G4ThreeVector the3Neutron = theNeutron.GetMomentum();
79  G4double nEnergy = theNeutron.GetTotalEnergy();
80  G4ThreeVector the3Target = theTarget.GetMomentum();
81  G4double tEnergy = theTarget.GetTotalEnergy();
82  G4ReactionProduct theCMS;
83  G4double totE = nEnergy+tEnergy;
84  G4ThreeVector the3CMS = the3Target+the3Neutron;
85  theCMS.SetMomentum(the3CMS);
86  G4double cmsMom = std::sqrt(the3CMS*the3CMS);
87  G4double sqrts = std::sqrt((totE-cmsMom)*(totE+cmsMom));
88  theCMS.SetMass(sqrts);
89  theCMS.SetTotalEnergy(totE);
90 
91  theNeutron.Lorentz(theNeutron, theCMS);
92  theTarget.Lorentz(theTarget, theCMS);
93  G4double en = theNeutron.GetTotalMomentum(); // already in CMS.
94  G4ThreeVector cms3Mom=theNeutron.GetMomentum(); // for neutron direction in CMS
95  G4double cms_theta=cms3Mom.theta();
96  G4double cms_phi=cms3Mom.phi();
97  G4ThreeVector tempVector;
98  tempVector.setX( std::cos(theta)*std::sin(cms_theta)*std::cos(cms_phi)
99  +std::sin(theta)*std::cos(phi)*std::cos(cms_theta)*std::cos(cms_phi)
100  -std::sin(theta)*std::sin(phi)*std::sin(cms_phi) );
101  tempVector.setY( std::cos(theta)*std::sin(cms_theta)*std::sin(cms_phi)
102  +std::sin(theta)*std::cos(phi)*std::cos(cms_theta)*std::sin(cms_phi)
103  +std::sin(theta)*std::sin(phi)*std::cos(cms_phi) );
104  tempVector.setZ( std::cos(theta)*std::cos(cms_theta)
105  -std::sin(theta)*std::cos(phi)*std::sin(cms_theta) );
106  tempVector *= en;
107  theNeutron.SetMomentum(tempVector);
108  theTarget.SetMomentum(-tempVector);
109  G4double tP = theTarget.GetTotalMomentum();
110  G4double tM = theTarget.GetMass();
111  theTarget.SetTotalEnergy(std::sqrt((tP+tM)*(tP+tM)-2.*tP*tM));
112 
113 
114  theNeutron.Lorentz(theNeutron, -1.*theCMS);
115 
116 //110913 Add Protection for very low energy (1e-6eV) scattering
117  if ( theNeutron.GetKineticEnergy() <= 0 )
118  {
119  theNeutron.SetTotalEnergy ( theNeutron.GetMass() * ( 1 + std::pow( 10 , -15.65 ) ) );
120  }
121 
122  theTarget.Lorentz(theTarget, -1.*theCMS);
123  if ( theTarget.GetKineticEnergy() < 0 )
124  {
125  theTarget.SetTotalEnergy ( theTarget.GetMass() * ( 1 + std::pow( 10 , -15.65 ) ) );
126  }
127 //110913 END
128 
129  theTarget.Lorentz(theTarget, -1.*theCMS);
130 
131  theResult->SetEnergyChange(theNeutron.GetKineticEnergy());
132  theResult->SetMomentumChange(theNeutron.GetMomentum().unit());
133  G4DynamicParticle* theRecoil = new G4DynamicParticle;
134 
135 // theRecoil->SetDefinition( ionTable->GetIon( iZ , iA ) );
136  theRecoil->SetDefinition( G4ParticleTable::GetParticleTable()->FindIon( iZ, iA , 0, iZ ));
137  theRecoil->SetMomentum( theTarget.GetMomentum() );
138 
139  theResult->AddSecondary( theRecoil );
140 
141  return theResult;
142 
143 }
144 
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
ThreeVector shoot(const G4int Ap, const G4int Af)
void SetMomentum(const G4ThreeVector &momentum)
G4double GetTotalMomentum() const
void Lorentz(const G4ReactionProduct &p1, const G4ReactionProduct &p2)
void SetKineticEnergy(const G4double en)
void SetMomentum(const G4double x, const G4double y, const G4double z)
int G4int
Definition: G4Types.hh:78
void setY(double)
void setZ(double)
void setX(double)
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
G4int GetNucleusEncoding(G4int iZ, G4int iA, G4int iM)
void SetMass(const G4double mas)
Hep3Vector vect() const
#define G4UniformRand()
Definition: Randomize.hh:87
float k_Boltzmann
Definition: hepunit.py:299
const G4ParticleDefinition * GetDefinition() const
double getElasticFinalState(double e_in, double temperature, double(*rng)(void *), void *rngState)
G4int Getm() const
Definition: G4Isotope.hh:100
std::map< G4int, G4LENDUsedTarget * > usedTarget_map
Definition: G4LENDModel.hh:79
G4double GetKineticEnergy() const
void SetTotalEnergy(const G4double en)
G4ErrorTarget * theTarget
Definition: errprop.cc:59
G4ParticleDefinition * FindIon(G4int atomicNumber, G4int atomicMass, G4double excitationEnergy)
double phi() const
const G4LorentzVector & Get4Momentum() const
double theta() const
G4double GetKineticEnergy() const
void SetEnergyChange(G4double anEnergy)
G4double GetTotalEnergy() const
G4double GetPDGMass() const
static G4ParticleTable * GetParticleTable()
const G4Isotope * GetIsotope()
Definition: G4Nucleus.hh:119
Hep3Vector unit() const
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
G4ThreeVector GetMomentum() const
G4double GetTemperature() const
Definition: G4Material.hh:180
const G4Material * GetMaterial() const
double G4double
Definition: G4Types.hh:76
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
G4LENDManager * lend_manager
Definition: G4LENDModel.hh:78
void SetMomentumChange(const G4ThreeVector &aV)
G4double GetMass() const
void AddSecondary(G4DynamicParticle *aP)