Geant4  10.01.p03
G4ParticleHPEnAngCorrelation.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 // particle_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 //
32 // P. Arce, June-2014 Conversion neutron_hp to particle_hp
33 //
35 #include "G4LorentzRotation.hh"
36 #include "G4LorentzVector.hh"
37 #include "G4RotationMatrix.hh"
38 #include "G4IonTable.hh"
39 
41 {
42  G4ReactionProduct * result = new G4ReactionProduct;
43 
44  // do we have an appropriate distribution
45  if(nProducts!=1) throw G4HadronicException(__FILE__, __LINE__, "More than one product in SampleOne");
46 
47  // get the result
48  G4ReactionProductVector * temp=0;
49  G4int i=0;
50  while(temp == 0) temp = theProducts[i++].Sample(anEnergy,1);
51 
52  // is the multiplicity correct
53  if(temp->size()!=1) throw G4HadronicException(__FILE__, __LINE__, "SampleOne: Yield not correct");
54 
55  // fill result
56  result = temp->operator[](0);
57 
58  // some garbage collection
59  delete temp;
60 
61  // return result
62  return result;
63 }
64 
66 {
68  G4int i;
70  G4ReactionProduct theCMS;
72 
73  if(frameFlag==2
74  || frameFlag==3) // Added for particle HP
75  {
76  // simplify and double check @
77  G4ThreeVector the3IncidentPart = theProjectileRP.GetMomentum(); //theProjectileRP has value in LAB
79  G4ThreeVector the3Target = theTarget.GetMomentum(); //theTarget has value in LAB
80  G4double tEnergy = theTarget.GetTotalEnergy();
81  G4double totE = nEnergy+tEnergy;
82  G4ThreeVector the3CMS = the3Target+the3IncidentPart;
83  theCMS.SetMomentum(the3CMS);
84  G4double cmsMom = std::sqrt(the3CMS*the3CMS);
85  G4double sqrts = std::sqrt((totE-cmsMom)*(totE+cmsMom));
86  theCMS.SetMass(sqrts);
87  theCMS.SetTotalEnergy(totE);
88  G4ReactionProduct aIncidentPart;
89  aIncidentPart.Lorentz(theProjectileRP, theCMS);
90  //TKDB 100413
91  //ENDF-6 Formats Manual ENDF-102
92  //CHAPTER 6. FILE 6: PRODUCT ENERGY-ANGLE DISTRIBUTIONS
93  //LCT Reference system for secondary energy and angle (incident energy is always given in the LAB system)
94  //anEnergy = aIncidentPart.GetKineticEnergy();
95  anEnergy = theProjectileRP.GetKineticEnergy(); //should be same argumment of "anEnergy"
96 
97  G4LorentzVector Ptmp (aIncidentPart.GetMomentum(), aIncidentPart.GetTotalEnergy());
98 
99  toZ.rotateZ(-1*Ptmp.phi());
100  toZ.rotateY(-1*Ptmp.theta());
101  }
103  G4LorentzRotation toLab(toZ.inverse()); //toLab only change axis NOT to LAB system
104  //- get first number of particles, to check if sum of Z and N is not bigger than target values
105  std::vector<int> nParticles;
106  bool bNPOK = true;
107  do {
108  G4int sumZ = 0;
109  G4int sumA = 0;
110  for(i=0; i<nProducts; i++)
111  {
112  G4int massCode = G4int(theProducts[i].GetMassCode());
113  G4int nPart;
114  nPart = theProducts[i].GetMultiplicity(anEnergy);
115  sumZ += massCode/1000 * nPart;
116  sumA += massCode % 1000 * nPart;
117 #ifdef G4PHPDEBUG
118  if( getenv("G4ParticleHPDebug") ) G4cout << i << " G4ParticleHPEnAngCorrelation::MULTIPLICITY " << massCode << " sumZ " << sumZ << " sumA " << sumA << " NPART " << nPart << G4endl;
119 #endif
120  nParticles.push_back( nPart );
121  }
122  bNPOK = true;
123  double targetZ = theTarget.GetDefinition()->GetAtomicNumber();
124  double targetA = theTarget.GetDefinition()->GetAtomicMass();
127  if( bAdjustFinalState ) {
128  if ( (sumZ != targetZ || sumA != targetA ) &&
129  (sumZ > targetZ || sumA > targetA
130  || ! G4IonTable::GetIonTable()->GetIon ( int(targetZ - sumZ), (int)(targetA - sumA), 0.0 ) ) ){ // e.g. Z=3, A=2
131  bNPOK = false;
132  nParticles.clear();
133 #ifdef G4PHPDEBUG
134  if( getenv("G4ParticleHPDebug") )
135  G4cerr << " WRONG MULTIPLICITY Z= " << sumZ
136  << " > " << targetZ
137  << " A= " << sumA
138  << " > " << targetA << G4endl;
139 #endif
140  }
141  }
142 
143  }while(!bNPOK);
144 
145  for(i=0; i<nProducts; i++)
146  {
147  //- if( nParticles[i] == 0 ) continue;
148  it = theProducts[i].Sample(anEnergy,nParticles[i]);
150  // if( getenv("G4PHPTEST") ) G4cout << " EnAnG energy sampled " << it->operator[](0)->GetKineticEnergy() << " aMeanEnergy " << aMeanEnergy << G4endl; // GDEB
151  if(aMeanEnergy>0)
152  {
153  theTotalMeanEnergy += aMeanEnergy;
154  }
155  else
156  {
157  theTotalMeanEnergy = anEnergy/nProducts+theProducts[i].GetQValue();
158  }
159  if(it!=0)
160  {
161  for(unsigned int ii=0; ii<it->size(); ii++)
162  {
163  //if(!getenv("G4PHP_NO_LORENTZ_BOOST")) {
164  G4LorentzVector pTmp1 (it->operator[](ii)->GetMomentum(),
165  it->operator[](ii)->GetTotalEnergy());
166  pTmp1 = toLab*pTmp1;
167  if( getenv("G4PHPTEST") ) G4cout << " G4particleHPEnAngCorrelation COS THETA " << std::cos(it->operator[](ii)->GetMomentum().theta()) << G4endl;
168  it->operator[](ii)->SetMomentum(pTmp1.vect());
169  it->operator[](ii)->SetTotalEnergy(pTmp1.e());
170  if( getenv("G4PHPTEST") ) G4cout << " G4particleHPEnAngCorrelation COS THETA after toLab " << std::cos(it->operator[](ii)->GetMomentum().theta()) << G4endl;
171 
172  if(frameFlag==1) // target rest //TK 100413 should be LAB?
173  {
174  it->operator[](ii)->Lorentz(*(it->operator[](ii)), -1.*theTarget); //TK 100413 Is this really need?
175  }
176  else if(frameFlag==2 ) // CMS
177  {
178 #ifdef G4PHPDEBUG
179  if( getenv("G4ParticleHPDebug") )
180  G4cout <<"G4ParticleHPEnAngCorrelation: before Lorentz boost "<<
181  it->at(ii)->GetKineticEnergy()<<" "<<
182  it->at(ii)->GetMomentum()<<G4endl;
183 #endif
184  it->operator[](ii)->Lorentz(*(it->operator[](ii)), -1.*theCMS);
185 #ifdef G4PHPDEBUG
186  if( getenv("G4ParticleHPDebug") )
187  G4cout <<"G4ParticleHPEnAngCorrelation: after Lorentz boost "<<
188  it->at(ii)->GetKineticEnergy()<<" "<<
189  it->at(ii)->GetMomentum()<<G4endl;
190 #endif
191  }
192  //TK120515 migrate frameFlag (MF6 LCT) = 3
193  else if(frameFlag==3) // CMS A<=4 other LAB
194  {
195  if ( theProducts[i].GetMassCode() > 4 ) //Alpha AWP 3.96713
196  {
197  //LAB
198  it->operator[](ii)->Lorentz(*(it->operator[](ii)), -1.*theTarget); //TK 100413 Is this really need?
199 #ifdef G4PHPDEBUG
200  if( getenv("G4ParticleHPDebug") )
201  G4cout <<"G4ParticleHPEnAngCorrelation: after Lorentz boost "<<
202  it->at(ii)->GetKineticEnergy()<<" "<<
203  it->at(ii)->GetMomentum()<<G4endl;
204 #endif
205  }
206  else
207  {
208  //CMS
209  it->operator[](ii)->Lorentz(*(it->operator[](ii)), -1.*theCMS);
210 #ifdef G4PHPDEBUG
211  if( getenv("G4ParticleHPDebug") )
212  G4cout <<"G4ParticleHPEnAngCorrelation: after Lorentz boost "<<
213  it->at(ii)->GetKineticEnergy()<<" "<<
214  it->at(ii)->GetMomentum()<<G4endl;
215 #endif
216  }
217  }
218  else
219  {
220  throw G4HadronicException(__FILE__, __LINE__, "G4ParticleHPEnAngCorrelation::Sample: The frame of the finalstate is not specified");
221  }
222  if( getenv("G4PHPTEST") ) G4cout << frameFlag << " G4particleHPEnAngCorrelation COS THETA after Lorentz " << std::cos(it->operator[](ii)->GetMomentum().theta()) << G4endl;
223 
224  // }//getenv("G4PHP_NO_LORENTZ_BOOST"))
225  // G4cout << ii << " EnAnG energy after boost " << it->operator[](ii)->GetKineticEnergy() << G4endl; //GDEB
226  result->push_back(it->operator[](ii));
227  }
228  delete it;
229  }
230  }
231 
232  return result;
233 }
234 
void Lorentz(const G4ReactionProduct &p1, const G4ReactionProduct &p2)
CLHEP::Hep3Vector G4ThreeVector
CLHEP::HepLorentzRotation G4LorentzRotation
void SetMomentum(const G4double x, const G4double y, const G4double z)
G4double MeanEnergyOfThisInteraction()
int G4int
Definition: G4Types.hh:78
G4ReactionProductVector * Sample(G4double anEnergy)
G4int GetAtomicNumber() const
std::vector< G4ReactionProduct * > G4ReactionProductVector
const G4ParticleDefinition * GetDefinition() const
void SetMass(const G4double mas)
G4GLOB_DLL std::ostream G4cout
void SetTotalEnergy(const G4double en)
G4int GetAtomicMass() const
static G4IonTable * GetIonTable()
Definition: G4IonTable.hh:78
G4double GetKineticEnergy() const
G4double GetTotalEnergy() const
G4int GetMultiplicity(G4double anEnergy)
G4ReactionProduct * SampleOne(G4double anEnergy)
G4ReactionProductVector * Sample(G4double anEnergy, G4int nParticles)
G4ThreeVector GetMomentum() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4GLOB_DLL std::ostream G4cerr
CLHEP::HepLorentzVector G4LorentzVector