Geant4  9.6.p02
 All Classes 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 
44  G4double ke = aTrack.GetKineticEnergy();
45 
46  //G4HadFinalState* theResult = new G4HadFinalState();
47  G4HadFinalState* theResult = &theParticleChange;
48  theResult->Clear();
49 
50  G4GIDI_target* aTarget = usedTarget_map.find( lend_manager->GetNucleusEncoding( iZ , iA ) )->second->GetTarget();
51  G4double aMu = aTarget->getElasticFinalState( ke*MeV, temp, NULL, NULL );
52 
54  G4double theta = std::acos( aMu );
55  //G4double sinth = std::sin( theta );
56 
57  G4ReactionProduct theNeutron( const_cast<G4ParticleDefinition *>( aTrack.GetDefinition() ) );
58  theNeutron.SetMomentum( aTrack.Get4Momentum().vect() );
59  theNeutron.SetKineticEnergy( ke );
60 
61 //G4cout << "iZ " << iZ << " iA " << iA << G4endl;
62 
63  G4ReactionProduct theTarget( G4ParticleTable::GetParticleTable()->FindIon( iZ , iA , 0 , iZ ) );
64 
66 
67 // add Thermal motion
68  G4double kT = k_Boltzmann*temp;
69  G4ThreeVector v ( G4RandGauss::shoot() * std::sqrt( kT*mass )
70  , G4RandGauss::shoot() * std::sqrt( kT*mass )
71  , G4RandGauss::shoot() * std::sqrt( kT*mass ) );
72  theTarget.SetMomentum( v );
73 
74  G4ThreeVector the3Neutron = theNeutron.GetMomentum();
75  G4double nEnergy = theNeutron.GetTotalEnergy();
76  G4ThreeVector the3Target = theTarget.GetMomentum();
77  G4double tEnergy = theTarget.GetTotalEnergy();
78  G4ReactionProduct theCMS;
79  G4double totE = nEnergy+tEnergy;
80  G4ThreeVector the3CMS = the3Target+the3Neutron;
81  theCMS.SetMomentum(the3CMS);
82  G4double cmsMom = std::sqrt(the3CMS*the3CMS);
83  G4double sqrts = std::sqrt((totE-cmsMom)*(totE+cmsMom));
84  theCMS.SetMass(sqrts);
85  theCMS.SetTotalEnergy(totE);
86 
87  theNeutron.Lorentz(theNeutron, theCMS);
88  theTarget.Lorentz(theTarget, theCMS);
89  G4double en = theNeutron.GetTotalMomentum(); // already in CMS.
90  G4ThreeVector cms3Mom=theNeutron.GetMomentum(); // for neutron direction in CMS
91  G4double cms_theta=cms3Mom.theta();
92  G4double cms_phi=cms3Mom.phi();
93  G4ThreeVector tempVector;
94  tempVector.setX( std::cos(theta)*std::sin(cms_theta)*std::cos(cms_phi)
95  +std::sin(theta)*std::cos(phi)*std::cos(cms_theta)*std::cos(cms_phi)
96  -std::sin(theta)*std::sin(phi)*std::sin(cms_phi) );
97  tempVector.setY( std::cos(theta)*std::sin(cms_theta)*std::sin(cms_phi)
98  +std::sin(theta)*std::cos(phi)*std::cos(cms_theta)*std::sin(cms_phi)
99  +std::sin(theta)*std::sin(phi)*std::cos(cms_phi) );
100  tempVector.setZ( std::cos(theta)*std::cos(cms_theta)
101  -std::sin(theta)*std::cos(phi)*std::sin(cms_theta) );
102  tempVector *= en;
103  theNeutron.SetMomentum(tempVector);
104  theTarget.SetMomentum(-tempVector);
105  G4double tP = theTarget.GetTotalMomentum();
106  G4double tM = theTarget.GetMass();
107  theTarget.SetTotalEnergy(std::sqrt((tP+tM)*(tP+tM)-2.*tP*tM));
108 
109 
110  theNeutron.Lorentz(theNeutron, -1.*theCMS);
111 
112 //110913 Add Protection for very low energy (1e-6eV) scattering
113  if ( theNeutron.GetKineticEnergy() <= 0 )
114  {
115  theNeutron.SetTotalEnergy ( theNeutron.GetMass() * ( 1 + std::pow( 10 , -15.65 ) ) );
116  }
117 
118  theTarget.Lorentz(theTarget, -1.*theCMS);
119  if ( theTarget.GetKineticEnergy() < 0 )
120  {
121  theTarget.SetTotalEnergy ( theTarget.GetMass() * ( 1 + std::pow( 10 , -15.65 ) ) );
122  }
123 //110913 END
124 
125  theTarget.Lorentz(theTarget, -1.*theCMS);
126 
127  theResult->SetEnergyChange(theNeutron.GetKineticEnergy());
128  theResult->SetMomentumChange(theNeutron.GetMomentum().unit());
129  G4DynamicParticle* theRecoil = new G4DynamicParticle;
130 
131 // theRecoil->SetDefinition( ionTable->GetIon( iZ , iA ) );
132  theRecoil->SetDefinition( G4ParticleTable::GetParticleTable()->FindIon( iZ, iA , 0, iZ ));
133  theRecoil->SetMomentum( theTarget.GetMomentum() );
134 
135  theResult->AddSecondary( theRecoil );
136 
137  return theResult;
138 
139 }
140