Geant4  10.01.p03
G4ParticleHPFinalState.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 //080721 Create adjust_final_state method by T. Koi
28 //080801 Residual reconstruction with theNDLDataA,Z (A, Z, and momentum are adjusted) by T. Koi
29 //101110 Set lower limit for gamma energy(1keV) by T. Koi
30 // P. Arce, June-2014 Conversion neutron_hp to particle_hp
31 //
32 
34 
35 #include "G4PhysicalConstants.hh"
36 #include "G4SystemOfUnits.hh"
37 #include "G4IonTable.hh"
38 #include "G4Gamma.hh"
39 #include "G4Neutron.hh"
40 #include "G4Proton.hh"
41 #include "G4Deuteron.hh"
42 #include "G4Triton.hh"
43 #include "G4Alpha.hh"
44 
45 
47 {
48 
49  G4double minimum_energy = 1*keV;
50 
51  if ( adjustResult != true ) return;
52 
53  G4int nSecondaries = theResult.GetNumberOfSecondaries();
54 
55  G4int sum_Z = 0;
56  G4int sum_A = 0;
57  G4int max_SecZ = 0;
58  G4int max_SecA = 0;
59  G4int imaxA = -1;
60  for ( int i = 0 ; i < nSecondaries ; i++ )
61  {
63  max_SecZ = std::max ( max_SecZ , theResult.GetSecondary( i )->GetParticle()->GetDefinition()->GetAtomicNumber() );
65  max_SecA = std::max ( max_SecA , theResult.GetSecondary( i )->GetParticle()->GetDefinition()->GetAtomicMass() );
66  if ( theResult.GetSecondary( i )->GetParticle()->GetDefinition()->GetAtomicMass() == max_SecA ) imaxA = i;
67 #ifdef G4PHPDEBUG
68  if( getenv("G4ParticleHPDebug")) G4cout << "G4ParticleHPFinalState::adjust_final_stat SECO " << i << " " <<theResult.GetSecondary( i )->GetParticle()->GetDefinition()->GetParticleName() << G4endl;
69 #endif
70 
71  }
72 
73  G4ParticleDefinition* resi_pd = NULL;
74 
75  G4double baseZNew = theBaseZ;
76  G4double baseANew = theBaseA;
78  baseANew ++;
79  } else if( theProjectile == G4Proton::Proton() ) {
80  baseZNew ++;
81  baseANew ++;
82  } else if( theProjectile == G4Deuteron::Deuteron() ) {
83  baseZNew ++;
84  baseANew += 2;
85  } else if( theProjectile == G4Triton::Triton() ) {
86  baseZNew ++;
87  baseANew += 3;
88  } else if( theProjectile == G4Alpha::Alpha() ) {
89  baseZNew += 2;
90  baseANew += 4;
91  }
92 
93 #ifdef G4PHPDEBUG
94  if( getenv("G4ParticleHPDebug")) G4cout << "G4ParticleHPFinalState::adjust_final_stat BaseZ " << baseZNew << " BaseA " << baseANew << " sum_Z " << sum_Z << " sum_A " << sum_A << G4endl;
95 #endif
96 
97  G4bool needOneMoreSec = false;
98  G4ParticleDefinition* oneMoreSec_pd = NULL;
99  if ( (int)(baseZNew - sum_Z) == 0 && (int)(baseANew - sum_A) == 0 )
100  {
101  //All secondaries are already created;
102  resi_pd = theResult.GetSecondary( imaxA )->GetParticle()->GetDefinition();
103  }
104  else
105  {
106  if ( max_SecA > int(baseANew - sum_A) )
107  {
108  //Most heavy secondary is interpreted as residual
109  resi_pd = theResult.GetSecondary( imaxA )->GetParticle()->GetDefinition();
110  needOneMoreSec = true;
111  }
112  else
113  {
114  //creation of residual is requierd
115  resi_pd = G4IonTable::GetIonTable()->GetIon ( int(baseZNew - sum_Z) , (int)(baseANew - sum_A) , 0.0 );
116  }
117 
118  if ( needOneMoreSec )
119  {
120  if ( int(baseZNew - sum_Z) == 0 && (int)(baseANew - sum_A) > 0 )
121  {
122  if ( int(baseANew - sum_A) > 1 ) G4cout << "More than one neutron is required for the balance of baryon number!" << G4endl;
123  oneMoreSec_pd = theProjectile;
124  }
125  else
126  {
127 #ifdef G4PHPDEBUG
128  if( getenv("G4ParticleHPDebug")) G4cout << this << "G4ParticleHPFinalState oneMoreSec_pd Z " << baseZNew << " - " << sum_Z << " A " << baseANew << " - " << sum_A << " projectile " << theProjectile->GetParticleName() << G4endl;
129 #endif
130  oneMoreSec_pd = G4IonTable::GetIonTable()->GetIon ( int(baseZNew - sum_Z) , (int)(baseANew - sum_A) , 0.0 );
131  if( !oneMoreSec_pd ) {
132  G4cerr << this << "G4ParticleHPFinalState oneMoreSec_pd Z " << baseZNew << " - " << sum_Z << " A " << baseANew << " - " << sum_A << " projectile " << theProjectile->GetParticleName() << G4endl;
133  G4Exception("G4ParticleHPFinalState:adjust_final_state",
134  "Warning",
135  JustWarning,
136  "No adjustment will be done!");
137  return;
138  }
139  }
140  }
141 
142  if ( resi_pd == NULL )
143  {
144  // theNDLDataZ,A has the Z and A of used NDL file
145  G4double ndlZNew = theNDLDataZ;
146  G4double ndlANew = theNDLDataA;
147  if( theProjectile == G4Neutron::Neutron() ) {
148  ndlANew ++;
149  } else if( theProjectile == G4Proton::Proton() ) {
150  ndlZNew ++;
151  ndlANew ++;
152  } else if( theProjectile == G4Deuteron::Deuteron() ) {
153  ndlZNew ++;
154  ndlANew += 2;
155  } else if( theProjectile == G4Triton::Triton() ) {
156  ndlZNew ++;
157  ndlANew += 3;
158  } else if( theProjectile == G4Alpha::Alpha() ) {
159  ndlZNew += 2;
160  ndlANew += 4;
161  }
162  // theNDLDataZ,A has the Z and A of used NDL file
163  if ( (int)(ndlZNew - sum_Z) == 0 && (int)(ndlANew - sum_A) == 0 )
164  {
165  G4int dif_Z = ( int ) ( theNDLDataZ - theBaseZ );
166  G4int dif_A = ( int ) ( theNDLDataA - theBaseA );
167  resi_pd = G4IonTable::GetIonTable()->GetIon ( max_SecZ - dif_Z , max_SecA - dif_A , 0.0 );
168  if( !resi_pd ) {
169  G4cerr << "G4ParticleHPFinalState resi_pd Z " << max_SecZ << " - " << dif_Z << " A " << max_SecA << " - " << dif_A << " projectile " << theProjectile->GetParticleName() << G4endl;
170  G4Exception("G4ParticleHPFinalState:adjust_final_state",
171  "Warning",
172  JustWarning,
173  "No adjustment will be done!");
174  return;
175  }
176 
177  for ( int i = 0 ; i < nSecondaries ; i++ )
178  {
179  if ( theResult.GetSecondary( i )->GetParticle()->GetDefinition()->GetAtomicNumber() == max_SecZ
180  && theResult.GetSecondary( i )->GetParticle()->GetDefinition()->GetAtomicMass() == max_SecA )
181  {
183  p = p * resi_pd->GetPDGMass()/ G4IonTable::GetIonTable()->GetIon ( max_SecZ , max_SecA , 0.0 )->GetPDGMass();
184  theResult.GetSecondary( i )->GetParticle()->SetDefinition( resi_pd );
186  }
187  }
188  }
189  }
190  }
191 
192 
193  G4LorentzVector secs_4p_lab( 0.0 );
194 
196  G4double fast = 0;
197  G4double slow = 1;
198  G4int ifast = 0;
199  G4int islow = 0;
200  G4int ires = -1;
201 
202  for ( G4int i = 0 ; i < n_sec ; i++ )
203  {
204 
205  //G4cout << "HP_DB " << i
206  // << " " << theResult.GetSecondary( i )->GetParticle()->GetDefinition()->GetParticleName()
207  // << " 4p " << theResult.GetSecondary( i )->GetParticle()->Get4Momentum()
208  // << " ke " << theResult.GetSecondary( i )->GetParticle()->Get4Momentum().e() - theResult.GetSecondary( i )->GetParticle()->GetDefinition()->GetPDGMass()
209  // << G4endl;
210 
211  secs_4p_lab += theResult.GetSecondary( i )->GetParticle()->Get4Momentum();
212 
213  G4double beta = 0;
215  {
216  beta = theResult.GetSecondary( i )->GetParticle()->Get4Momentum().beta();
217  }
218  else
219  {
220  beta = 1;
221  }
222 
223  if ( theResult.GetSecondary( i )->GetParticle()->GetDefinition() == resi_pd ) ires = i;
224 
225  if ( slow > beta && beta != 0 )
226  {
227  slow = beta;
228  islow = i;
229  }
230 
231  if ( fast <= beta )
232  {
233  if ( fast != 1 )
234  {
235  fast = beta;
236  ifast = i;
237  }
238  else
239  {
240 // fast is already photon then check E
242  if ( e > theResult.GetSecondary( ifast )->GetParticle()->Get4Momentum().e() )
243  {
244 // among photons, the highest E becomes the fastest
245  ifast = i;
246  }
247  }
248  }
249  }
250 
251 
252  G4LorentzVector dif_4p = init_4p_lab - secs_4p_lab;
253 
254  //G4cout << "HP_DB dif_4p " << init_4p_lab - secs_4p_lab << G4endl;
255  //G4cout << "HP_DB dif_3p mag " << ( dif_4p.v() ).mag() << G4endl;
256  //G4cout << "HP_DB dif_e " << dif_4p.e() - ( dif_4p.v() ).mag()<< G4endl;
257 
258  G4LorentzVector p4(0);
259  if ( ires == -1 )
260  {
261 // Create and Add Residual Nucleus
262  ires = nSecondaries;
263  nSecondaries += 1;
264 
265  G4DynamicParticle* res = new G4DynamicParticle ( resi_pd , dif_4p.v() );
266  theResult.AddSecondary ( res );
267 
268  p4 = res->Get4Momentum();
269  if ( slow > p4.beta() )
270  {
271  slow = p4.beta();
272  islow = ires;
273  }
274  dif_4p = init_4p_lab - ( secs_4p_lab + p4 );
275  }
276 
277  if ( needOneMoreSec && oneMoreSec_pd)
278  //
279  // fca: this is not a fix, this is a crash avoidance...
280  // fca: the baryon number is still wrong, most probably because it
281  // fca: should have been decreased, but since we could not create a particle
282  // fca: we just do not add it
283  //
284  {
285  nSecondaries += 1;
286  G4DynamicParticle* one = new G4DynamicParticle ( oneMoreSec_pd , dif_4p.v() );
287  theResult.AddSecondary ( one );
288  p4 = one->Get4Momentum();
289  if ( slow > p4.beta() )
290  {
291  slow = p4.beta();
292  islow = nSecondaries-1; //Because the first is 0th, so the last becomes "nSecondaries-1"
293  }
294  dif_4p = init_4p_lab - ( secs_4p_lab + p4 );
295  }
296 
297  //Which is bigger dif_p or dif_e
298 
299  if ( dif_4p.v().mag() < std::abs( dif_4p.e() ) )
300  {
301 
302  // Adjust p
303  //if ( dif_4p.v().mag() < 1*MeV )
304  if ( minimum_energy < dif_4p.v().mag() && dif_4p.v().mag() < 1*MeV )
305  {
306 
307  nSecondaries += 1;
308  theResult.AddSecondary ( new G4DynamicParticle ( G4Gamma::Gamma() , dif_4p.v() ) );
309 
310  }
311  else
312  {
313  //G4cout << "HP_DB Difference in dif_p is too large (>1MeV) or too small(<1keV) to adjust, so that give up tuning" << G4endl;
314  }
315 
316  }
317  else
318  {
319 
320  // dif_p > dif_e
321  // at first momentum
322  // Move residual momentum
323 
324  p4 = theResult.GetSecondary( ires )->GetParticle()->Get4Momentum();
325  theResult.GetSecondary( ires )->GetParticle()->SetMomentum( p4.v() + dif_4p.v() );
326  dif_4p = init_4p_lab - ( secs_4p_lab - p4 + theResult.GetSecondary( ires )->GetParticle()->Get4Momentum() );
327 
328  //G4cout << "HP_DB new residual kinetic energy " << theResult.GetSecondary( ires )->GetParticle()->GetKineticEnergy() << G4endl;
329 
330  }
331 
332  G4double dif_e = dif_4p.e() - ( dif_4p.v() ).mag();
333  //G4cout << "HP_DB dif_e " << dif_e << G4endl;
334 
335  if ( dif_e > 0 )
336  {
337 
338 // create 2 gamma
339 
340  nSecondaries += 2;
341  G4double e1 = ( dif_4p.e() -dif_4p.v().mag() ) / 2;
342 
343  if ( minimum_energy < e1 )
344  {
345  G4double costh = 2.*G4UniformRand()-1.;
346  G4double phi = twopi*G4UniformRand();
347  G4ThreeVector dir( std::sin(std::acos(costh))*std::cos(phi),
348  std::sin(std::acos(costh))*std::sin(phi),
349  costh);
350  theResult.AddSecondary ( new G4DynamicParticle ( G4Gamma::Gamma() , e1*dir ) );
351  theResult.AddSecondary ( new G4DynamicParticle ( G4Gamma::Gamma() , -e1*dir ) );
352  }
353  else
354  {
355  //G4cout << "HP_DB Difference is too small(<1keV) to adjust, so that neglect it" << G4endl;
356  }
357 
358  }
359  else //dif_e < 0
360  {
361 
362 // At first reduce KE of the fastest secondary;
365  G4ThreeVector dir = ( theResult.GetSecondary( ifast )->GetParticle()->GetMomentum() ).unit();
366 
367  //G4cout << "HP_DB ifast " << ifast << " ke0 " << ke0 << G4endl;
368 
369  if ( ke0 + dif_e > 0 )
370  {
371  theResult.GetSecondary( ifast )->GetParticle()->SetKineticEnergy( ke0 + dif_e );
372  G4ThreeVector dp = p0 - theResult.GetSecondary( ifast )->GetParticle()->GetMomentum();
373 
375  //theResult.GetSecondary( islow )->GetParticle()->SetMomentum( p - dif_e*dir );
376  theResult.GetSecondary( islow )->GetParticle()->SetMomentum( p + dp );
377  }
378  else
379  {
380  //G4cout << "HP_DB Difference in dif_e too large ( <0MeV ) to adjust, so that give up tuning" << G4endl;
381  }
382 
383  }
384 
385 }
386 
static const double MeV
Definition: G4SIunits.hh:193
void SetMomentum(const G4ThreeVector &momentum)
G4HadSecondary * GetSecondary(size_t i)
G4double GetKineticEnergy() const
CLHEP::Hep3Vector G4ThreeVector
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:463
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
G4int GetAtomicNumber() const
G4ParticleDefinition * theProjectile
void adjust_final_state(G4LorentzVector)
#define G4UniformRand()
Definition: Randomize.hh:93
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
static G4Triton * Triton()
Definition: G4Triton.cc:95
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:94
G4int GetAtomicMass() const
void SetKineticEnergy(G4double aEnergy)
G4LorentzVector Get4Momentum() const
static G4IonTable * GetIonTable()
Definition: G4IonTable.hh:78
static const G4double e1
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
const G4double p0
G4double GetPDGMass() const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4DynamicParticle * GetParticle()
#define G4endl
Definition: G4ios.hh:61
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
static const double keV
Definition: G4SIunits.hh:195
double G4double
Definition: G4Types.hh:76
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
G4int GetNumberOfSecondaries() const
G4ThreeVector GetMomentum() const
G4GLOB_DLL std::ostream G4cerr
CLHEP::HepLorentzVector G4LorentzVector