Geant4  10.02.p02
G4ParticleHPDiscreteTwoBody.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 //080612 Bug fix contribution from Benoit Pirard and Laurent Desorgher (Univ. Bern) #2,3
31 //080709 Bug fix Sampling Legendre expansion by T. Koi
32 //101110 Bug fix in MF=6, LAW=2 case; contribution from E. Mendoza, D. Cano-Ott (CIEMAT)
33 //
34 // P. Arce, June-2014 Conversion neutron_hp to particle_hp
35 //
37 #include "G4Gamma.hh"
38 #include "G4Electron.hh"
39 #include "G4Positron.hh"
40 #include "G4Neutron.hh"
41 #include "G4Proton.hh"
42 #include "G4Deuteron.hh"
43 #include "G4Triton.hh"
44 #include "G4He3.hh"
45 #include "G4Alpha.hh"
46 #include "G4ParticleHPVector.hh"
48 
50 { // Interpolation still only for the most used parts; rest to be Done @@@@@
51  G4ReactionProduct * result = new G4ReactionProduct;
52  G4int Z = static_cast<G4int>(massCode/1000);
53  G4int A = static_cast<G4int>(massCode-1000*Z);
54 
55  if(massCode==0)
56  {
57  result->SetDefinition(G4Gamma::Gamma());
58  }
59  else if(A==0)
60  {
62  if(Z==1) result->SetDefinition(G4Positron::Positron());
63  }
64  else if(A==1)
65  {
67  if(Z==1) result->SetDefinition(G4Proton::Proton());
68  }
69  else if(A==2)
70  {
72  }
73  else if(A==3)
74  {
75  result->SetDefinition(G4Triton::Triton());
76  if(Z==2) result->SetDefinition(G4He3::He3());
77  }
78  else if(A==4)
79  {
80  result->SetDefinition(G4Alpha::Alpha());
81  if(Z!=2) throw G4HadronicException(__FILE__, __LINE__, "Unknown ion case 1");
82  }
83  else
84  {
85  throw G4HadronicException(__FILE__, __LINE__, "G4ParticleHPDiscreteTwoBody: Unknown ion case 2");
86  }
87 
88 // get cosine(theta)
89  G4int i(0), it(0);
90  G4double cosTh(0);
91  for(i=0; i<nEnergy; i++)
92  {
93  it = i;
94  if(theCoeff[i].GetEnergy()>anEnergy) break;
95  }
96  if(it==0||it==nEnergy-1)
97  {
98  if(theCoeff[it].GetRepresentation()==0)
99  {
100 //TK Legendre expansion
101  G4ParticleHPLegendreStore theStore(1);
102  theStore.SetCoeff(0, theCoeff);
103  theStore.SetManager(theManager);
104  //cosTh = theStore.SampleMax(anEnergy);
105  //080612TK contribution from Benoit Pirard and Laurent Desorgher (Univ. Bern) #3
106  cosTh = theStore.SampleDiscreteTwoBody(anEnergy);
107  }
108  else if(theCoeff[it].GetRepresentation()==12) // means LINLIN
109  {
110  G4ParticleHPVector theStore;
111  G4InterpolationManager aManager;
112  aManager.Init(LINLIN, theCoeff[it].GetNumberOfPoly()/2);
113  theStore.SetInterpolationManager(aManager);
114  for(i=0;i<theCoeff[it].GetNumberOfPoly(); i++)
115  {
116  //101110
117  //theStore.SetX(i, theCoeff[it].GetCoeff(i));
118  //theStore.SetY(i, theCoeff[it].GetCoeff(i));
119  theStore.SetX(i/2, theCoeff[it].GetCoeff(i));
120  theStore.SetY(i/2, theCoeff[it].GetCoeff(i+1));
121  i++;
122  }
123  cosTh = theStore.Sample();
124  }
125  else if(theCoeff[it].GetRepresentation()==14) //this is LOGLIN
126  {
127  G4ParticleHPVector theStore;
128  G4InterpolationManager aManager;
129  aManager.Init(LOGLIN, theCoeff[it].GetNumberOfPoly()/2);
130  theStore.SetInterpolationManager(aManager);
131  for(i=0;i<theCoeff[it].GetNumberOfPoly(); i++)
132  {
133  //101110
134  //theStore.SetX(i, theCoeff[it].GetCoeff(i));
135  //theStore.SetY(i, theCoeff[it].GetCoeff(i));
136  theStore.SetX(i/2, theCoeff[it].GetCoeff(i));
137  theStore.SetY(i/2, theCoeff[it].GetCoeff(i+1));
138  i++;
139  }
140  cosTh = theStore.Sample();
141  }
142  else
143  {
144  throw G4HadronicException(__FILE__, __LINE__, "unknown representation type in Two-body scattering");
145  }
146  }
147  else
148  {
149  if(!bCheckDiffCoeffRepr || theCoeff[it].GetRepresentation() == theCoeff[it-1].GetRepresentation())
150  {
151  if(theCoeff[it].GetRepresentation()==0)
152  {
153 //TK Legendre expansion
154  G4ParticleHPLegendreStore theStore(2);
155  theStore.SetCoeff(0, &(theCoeff[it-1]));
156  theStore.SetCoeff(1, &(theCoeff[it]));
157  G4InterpolationManager aManager;
158  aManager.Init(theManager.GetScheme(it), 2);
159  theStore.SetManager(aManager);
160  //cosTh = theStore.SampleMax(anEnergy);
161 //080709 TKDB
162  cosTh = theStore.SampleDiscreteTwoBody(anEnergy);
163  }
164  else if(theCoeff[it].GetRepresentation()==12) // LINLIN
165  {
166  G4ParticleHPVector theBuff1;
167  G4InterpolationManager aManager1;
168  aManager1.Init(LINLIN, theCoeff[it-1].GetNumberOfPoly()/2);
169  theBuff1.SetInterpolationManager(aManager1);
170  for(i=0;i<theCoeff[it-1].GetNumberOfPoly(); i++)
171  {
172  //101110
173  //theBuff1.SetX(i, theCoeff[it-1].GetCoeff(i));
174  //theBuff1.SetY(i, theCoeff[it-1].GetCoeff(i));
175  theBuff1.SetX(i/2, theCoeff[it-1].GetCoeff(i));
176  theBuff1.SetY(i/2, theCoeff[it-1].GetCoeff(i+1));
177  i++;
178  }
179  G4ParticleHPVector theBuff2;
180  G4InterpolationManager aManager2;
181  aManager2.Init(LINLIN, theCoeff[it].GetNumberOfPoly()/2);
182  theBuff2.SetInterpolationManager(aManager2);
183  for(i=0;i<theCoeff[it].GetNumberOfPoly(); i++)
184  {
185  //theBuff2.SetX(i, theCoeff[it].GetCoeff(i));
186  //theBuff2.SetY(i, theCoeff[it].GetCoeff(i));
187  theBuff2.SetX(i/2, theCoeff[it].GetCoeff(i));
188  theBuff2.SetY(i/2, theCoeff[it].GetCoeff(i+1));
189  i++;
190  }
191 
192  G4double x1 = theCoeff[it-1].GetEnergy();
193  G4double x2 = theCoeff[it].GetEnergy();
194  G4double x = anEnergy;
195  G4double y1, y2, y, mu;
196 
197  G4ParticleHPVector theStore1;
198  theStore1.SetInterpolationManager(aManager1);
199  G4ParticleHPVector theStore2;
200  theStore2.SetInterpolationManager(aManager2);
201  G4ParticleHPVector theStore;
202 
203  // for fixed mu get p1, p2 and interpolate according to x
204  for(i=0; i<theBuff1.GetVectorLength(); i++)
205  {
206  mu = theBuff1.GetX(i);
207  y1 = theBuff1.GetY(i);
208  y2 = theBuff2.GetY(mu);
209  y = theInt.Interpolate(theManager.GetScheme(it), x, x1,x2,y1,y2);
210  theStore1.SetData(i, mu, y);
211  }
212  for(i=0; i<theBuff2.GetVectorLength(); i++)
213  {
214  mu = theBuff2.GetX(i);
215  y1 = theBuff2.GetY(i);
216  y2 = theBuff1.GetY(mu);
217  y = theInt.Interpolate(theManager.GetScheme(it), x, x1,x2,y1,y2);
218  theStore2.SetData(i, mu, y);
219  }
220  theStore.Merge(&theStore1, &theStore2); // merge takes care of interpolationschemes
221  cosTh = theStore.Sample();
222  }
223  else if(theCoeff[it].GetRepresentation()==14) //TK LOG_LIN
224  {
225  G4ParticleHPVector theBuff1;
226  G4InterpolationManager aManager1;
227  aManager1.Init(LOGLIN, theCoeff[it-1].GetNumberOfPoly()/2);
228  theBuff1.SetInterpolationManager(aManager1);
229  for(i=0;i<theCoeff[it-1].GetNumberOfPoly(); i++)
230  {
231  //101110
232  //theBuff1.SetX(i, theCoeff[it-1].GetCoeff(i));
233  //theBuff1.SetY(i, theCoeff[it-1].GetCoeff(i));
234  theBuff1.SetX(i/2, theCoeff[it-1].GetCoeff(i));
235  theBuff1.SetY(i/2, theCoeff[it-1].GetCoeff(i+1));
236  i++;
237  }
238 
239  G4ParticleHPVector theBuff2;
240  G4InterpolationManager aManager2;
241  aManager2.Init(LOGLIN, theCoeff[it].GetNumberOfPoly()/2);
242  theBuff2.SetInterpolationManager(aManager2);
243  for(i=0;i<theCoeff[it].GetNumberOfPoly(); i++)
244  {
245  //101110
246  //theBuff2.SetX(i, theCoeff[it].GetCoeff(i));
247  //theBuff2.SetY(i, theCoeff[it].GetCoeff(i));
248  theBuff2.SetX(i/2, theCoeff[it].GetCoeff(i));
249  theBuff2.SetY(i/2, theCoeff[it].GetCoeff(i+1));
250  i++;
251  }
252 
253  G4double x1 = theCoeff[it-1].GetEnergy();
254  G4double x2 = theCoeff[it].GetEnergy();
255  G4double x = anEnergy;
256  G4double y1, y2, y, mu;
257 
258  G4ParticleHPVector theStore1;
259  theStore1.SetInterpolationManager(aManager1);
260  G4ParticleHPVector theStore2;
261  theStore2.SetInterpolationManager(aManager2);
262  G4ParticleHPVector theStore;
263 
264  // for fixed mu get p1, p2 and interpolate according to x
265  for(i=0; i<theBuff1.GetVectorLength(); i++)
266  {
267  mu = theBuff1.GetX(i);
268  y1 = theBuff1.GetY(i);
269  y2 = theBuff2.GetY(mu);
270  y = theInt.Interpolate(theManager.GetScheme(it), x, x1,x2,y1,y2);
271  theStore1.SetData(i, mu, y);
272  }
273  for(i=0; i<theBuff2.GetVectorLength(); i++)
274  {
275  mu = theBuff2.GetX(i);
276  y1 = theBuff2.GetY(i);
277  y2 = theBuff1.GetY(mu);
278  y = theInt.Interpolate(theManager.GetScheme(it), x, x1,x2,y1,y2);
279  theStore2.SetData(i, mu, y);
280  }
281  theStore.Merge(&theStore1, &theStore2);
282  cosTh = theStore.Sample();
283  }
284  else
285  {
286  throw G4HadronicException(__FILE__, __LINE__, "Two neighbouring distributions with different interpolation");
287  }
288  }
289  else
290  {
291  G4cout << " theCoeff[it].GetRepresent MEM " << it << " " << &theCoeff[it] << " " << &theCoeff[it-1] << G4endl;
292  G4cout << " theCoeff[it].GetRepresent " << it << " " << theCoeff[it].GetRepresentation() << " != " << theCoeff[it-1].GetRepresentation() << G4endl;
293 
294  throw G4HadronicException(__FILE__, __LINE__, "unknown representation type in Two-body scattering, case 2");
295  }
296  }
297 
298 // now get the energy from kinematics and Q-value.
299 
300  //G4double restEnergy = anEnergy+GetQValue();
301 
302 // assumed to be in CMS @@@@@@@@@@@@@@@@@
303 
304  //080612TK contribution from Benoit Pirard and Laurent Desorgher (Univ. Bern) #2
305  //G4double residualMass = GetTarget()->GetMass() + GetNeutron()->GetMass()
306  // - result->GetMass() - GetQValue();
307  //G4double kinE = restEnergy/(1+result->GetMass()/residualMass); // non relativistic @@
309  G4double A1prim = result->GetMass()/GetProjectileRP()->GetMass();
310  //G4double E1 = (A1+1)*(A1+1)/A1/A1*anEnergy;
311  //Bug fix Bugzilla #1815
312  G4double E1 = anEnergy;
313  G4double kinE = (A1+1-A1prim)/(A1+1)/(A1+1)*(A1*E1+(1+A1)*GetQValue());
314 
315  result->SetKineticEnergy(kinE); // non relativistic @@
317  G4double theta = std::acos(cosTh);
318  G4double sinth = std::sin(theta);
319  G4double mtot = result->GetTotalMomentum();
320  G4ThreeVector tempVector(mtot*sinth*std::cos(phi), mtot*sinth*std::sin(phi), mtot*std::cos(theta) );
321  result->SetMomentum(tempVector);
322 
323 // some garbage collection
324 
325 // return the result
326  return result;
327 }
G4int GetVectorLength() const
G4double GetTotalMomentum() const
CLHEP::Hep3Vector G4ThreeVector
void Init(G4int aScheme, G4int aRange)
void SetKineticEnergy(const G4double en)
void SetMomentum(const G4double x, const G4double y, const G4double z)
void SetData(G4int i, G4double x, G4double y)
G4double SampleDiscreteTwoBody(G4double anEnergy)
void SetInterpolationManager(const G4InterpolationManager &aManager)
int G4int
Definition: G4Types.hh:78
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
void SetY(G4int i, G4double x)
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
double A(double temperature)
G4InterpolationScheme GetScheme(G4int index) const
static const double twopi
Definition: G4SIunits.hh:75
static G4Triton * Triton()
Definition: G4Triton.cc:95
void Merge(G4ParticleHPVector *active, G4ParticleHPVector *passive)
static G4Proton * Proton()
Definition: G4Proton.cc:93
G4double GetX(G4int i) const
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
G4double GetY(G4double x)
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:94
G4double Interpolate(G4InterpolationScheme aScheme, G4double x, G4double x1, G4double x2, G4double y1, G4double y2) const
void SetManager(G4InterpolationManager &aManager)
static G4Positron * Positron()
Definition: G4Positron.cc:94
void SetCoeff(G4int i, G4int l, G4double coeff)
const G4double x[NPOINTSGL]
G4ReactionProduct * Sample(G4double anEnergy, G4double massCode, G4double mass)
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
double G4double
Definition: G4Types.hh:76
static G4He3 * He3()
Definition: G4He3.cc:94
void SetX(G4int i, G4double e)
G4double GetMass() const
G4ParticleHPLegendreTable * theCoeff