Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4NeutronHPEnAngCorrelation.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// neutron_hp -- source file
27// J.P. Wellisch, Nov-1996
28// A prototype of the low energy neutron transport model.
29//
30// 100413 Fix bug in incidence energy by T. Koi
31//
33#include "G4LorentzRotation.hh"
34#include "G4LorentzVector.hh"
35#include "G4RotationMatrix.hh"
36
38{
40
41 // do we have an appropriate distribution
42 if(nProducts!=1) throw G4HadronicException(__FILE__, __LINE__, "More than one product in SampleOne");
43
44 // get the result
46 G4int i=0;
47 while(temp == 0) temp = theProducts[i++].Sample(anEnergy);
48
49 // is the multiplicity correct
50 if(temp->size()!=1) throw G4HadronicException(__FILE__, __LINE__, "SampleOne: Yield not correct");
51
52 // fill result
53 result = temp->operator[](0);
54
55 // some garbage collection
56 delete temp;
57
58 // return result
59 return result;
60}
61
63{
65 G4int i;
67 G4ReactionProduct theCMS;
69 //TK120515 migrate frameFlag (MF6 LCT) = 3
70 //if(frameFlag==2)
71 if(frameFlag==2||frameFlag==3)
72 {
73 // simplify and double check @
74 G4ThreeVector the3Neutron = theNeutron.GetMomentum(); //theNeutron has value in LAB
75 G4double nEnergy = theNeutron.GetTotalEnergy();
76 G4ThreeVector the3Target = theTarget.GetMomentum(); //theTarget has value in LAB
77 G4double tEnergy = theTarget.GetTotalEnergy();
78 G4double totE = nEnergy+tEnergy;
79 G4ThreeVector the3CMS = the3Target+the3Neutron;
80 theCMS.SetMomentum(the3CMS);
81 G4double cmsMom = std::sqrt(the3CMS*the3CMS);
82 G4double sqrts = std::sqrt((totE-cmsMom)*(totE+cmsMom));
83 theCMS.SetMass(sqrts);
84 theCMS.SetTotalEnergy(totE);
85 G4ReactionProduct aNeutron;
86 aNeutron.Lorentz(theNeutron, theCMS);
87 //TKDB 100413
88 //ENDF-6 Formats Manual ENDF-102
89 //CHAPTER 6. FILE 6: PRODUCT ENERGY-ANGLE DISTRIBUTIONS
90 //LCT Reference system for secondary energy and angle (incident energy is always given in the LAB system)
91 //anEnergy = aNeutron.GetKineticEnergy();
92 anEnergy = theNeutron.GetKineticEnergy(); //should be same argumment of "anEnergy"
93
94 G4LorentzVector Ptmp (aNeutron.GetMomentum(), aNeutron.GetTotalEnergy());
95
96 toZ.rotateZ(-1*Ptmp.phi());
97 toZ.rotateY(-1*Ptmp.theta());
98 }
99 theTotalMeanEnergy=0;
100 G4LorentzRotation toLab(toZ.inverse()); //toLab only change axis NOT to LAB system
101 for(i=0; i<nProducts; i++)
102 {
103 it = theProducts[i].Sample(anEnergy);
104 G4double aMeanEnergy = theProducts[i].MeanEnergyOfThisInteraction();
105 if(aMeanEnergy>0)
106 {
107 theTotalMeanEnergy += aMeanEnergy;
108 }
109 else
110 {
111 theTotalMeanEnergy = anEnergy/nProducts+theProducts[i].GetQValue();
112 }
113 if(it!=0)
114 {
115 for(unsigned int ii=0; ii<it->size(); ii++)
116 {
117 G4LorentzVector pTmp1 (it->operator[](ii)->GetMomentum(),
118 it->operator[](ii)->GetTotalEnergy());
119 pTmp1 = toLab*pTmp1;
120 it->operator[](ii)->SetMomentum(pTmp1.vect());
121 it->operator[](ii)->SetTotalEnergy(pTmp1.e());
122 if(frameFlag==1) // target rest //TK 100413 should be LAB?
123 {
124 it->operator[](ii)->Lorentz(*(it->operator[](ii)), -1.*theTarget); //TK 100413 Is this really need?
125 }
126 else if(frameFlag==2) // CMS
127 {
128#ifdef G4_NHP_DEBUG
129 cout <<"G4NeutronHPEnAngCorrelation: "<<
130 it->at(ii)->GetTotalEnergy()<<" "<<
131 it->at(ii)->GetMomentum()<<G4endl;
132#endif
133 it->operator[](ii)->Lorentz(*(it->operator[](ii)), -1.*theCMS);
134 }
135 //TK120515 migrate frameFlag (MF6 LCT) = 3
136 else if(frameFlag==3) // CMS A<=4 other LAB
137 {
138 if ( theProducts[i].GetMassCode() > 4 ) //Alpha AWP 3.96713
139 {
140 //LAB
141 it->operator[](ii)->Lorentz(*(it->operator[](ii)), -1.*theTarget); //TK 100413 Is this really need?
142 }
143 else
144 {
145 //CMS
146 it->operator[](ii)->Lorentz(*(it->operator[](ii)), -1.*theCMS);
147 }
148 }
149 else
150 {
151 throw G4HadronicException(__FILE__, __LINE__, "G4NeutronHPEnAngCorrelation::Sample: The frame of the finalstate is not specified");
152 }
153 result->push_back(it->operator[](ii));
154 }
155 delete it;
156 }
157 }
158 return result;
159}
160
std::vector< G4ReactionProduct * > G4ReactionProductVector
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
HepLorentzRotation & rotateY(double delta)
HepLorentzRotation & rotateZ(double delta)
HepLorentzRotation inverse() const
double theta() const
Hep3Vector vect() const
G4ReactionProduct * SampleOne(G4double anEnergy)
G4ReactionProductVector * Sample(G4double anEnergy)
G4double MeanEnergyOfThisInteraction()
G4ReactionProductVector * Sample(G4double anEnergy)
void SetMomentum(const G4double x, const G4double y, const G4double z)
void SetTotalEnergy(const G4double en)
G4double GetKineticEnergy() const
G4double GetTotalEnergy() const
G4ThreeVector GetMomentum() const
void Lorentz(const G4ReactionProduct &p1, const G4ReactionProduct &p2)
void SetMass(const G4double mas)