Geant4  10.01.p03
G4NeutronHPFinalState.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 
31 #include "G4NeutronHPFinalState.hh"
32 #include "G4PhysicalConstants.hh"
33 #include "G4SystemOfUnits.hh"
34 #include "G4IonTable.hh"
35 #include "G4Gamma.hh"
36 #include "G4Neutron.hh"
37 
38 //G4ThreadLocal G4HadFinalState G4NeutronHPFinalState::theResult;
39 
41 {
42 
43  G4double minimum_energy = 1*keV;
44 
45  //if ( adjustResult != true ) return;
46  if ( G4NeutronHPManager::GetInstance()->GetDoNotAdjustFinalState() ) return;
47 
48  G4int nSecondaries = theResult.Get()->GetNumberOfSecondaries();
49 
50  G4int sum_Z = 0;
51  G4int sum_A = 0;
52  G4int max_SecZ = 0;
53  G4int max_SecA = 0;
54  G4int imaxA = -1;
55  for ( int i = 0 ; i < nSecondaries ; i++ )
56  {
58  max_SecZ = std::max ( max_SecZ , theResult.Get()->GetSecondary( i )->GetParticle()->GetDefinition()->GetAtomicNumber() );
60  max_SecA = std::max ( max_SecA , theResult.Get()->GetSecondary( i )->GetParticle()->GetDefinition()->GetAtomicMass() );
61  if ( theResult.Get()->GetSecondary( i )->GetParticle()->GetDefinition()->GetAtomicMass() == max_SecA ) imaxA = i;
62  }
63 
64  G4ParticleDefinition* resi_pd = NULL;
65  G4bool needOneMoreSec = false;
66  G4ParticleDefinition* oneMoreSec_pd = NULL;
67  if ( (int)(theBaseZ - sum_Z) == 0 && (int)(theBaseA + 1 - sum_A) == 0 )
68  {
69  //All secondaries are already created;
70  resi_pd = theResult.Get()->GetSecondary( imaxA )->GetParticle()->GetDefinition();
71  }
72  else
73  {
74  if ( max_SecA > int(theBaseA + 1 - sum_A) )
75  {
76  //Most heavy secondary is interpreted as residual
77  resi_pd = theResult.Get()->GetSecondary( imaxA )->GetParticle()->GetDefinition();
78  needOneMoreSec = true;
79  }
80  else
81  {
82  //creation of residual is requierd
83  resi_pd = G4IonTable::GetIonTable()->GetIon ( int(theBaseZ - sum_Z) , (int)(theBaseA + 1 - sum_A) , 0.0 );
84  }
85 
86  if ( needOneMoreSec )
87  {
88  if ( int(theBaseZ - sum_Z) == 0 && (int)(theBaseA + 1 - sum_A) > 0 )
89  {
90  if ( int(theBaseA + 1 - sum_A) > 1 ) G4cout << "More than one neutron is required for the balance of baryon number!" << G4endl;
91  oneMoreSec_pd = G4Neutron::Neutron();
92  }
93  else
94  {
95  oneMoreSec_pd = G4IonTable::GetIonTable()->GetIon ( int(theBaseZ - sum_Z) , (int)(theBaseA + 1 - sum_A) , 0.0 );
96  }
97  }
98 
99  if ( resi_pd == NULL )
100  {
101  // theNDLDataZ,A has the Z and A of used NDL file
102  if ( (int)(theNDLDataZ - sum_Z) == 0 && (int)(theNDLDataA + 1 - sum_A) == 0 )
103  {
104  G4int dif_Z = ( int ) ( theNDLDataZ - theBaseZ );
105  G4int dif_A = ( int ) ( theNDLDataA - theBaseA );
106  resi_pd = G4IonTable::GetIonTable()->GetIon ( max_SecZ - dif_Z , max_SecA - dif_A , 0.0 );
107  for ( int i = 0 ; i < nSecondaries ; i++ )
108  {
109  if ( theResult.Get()->GetSecondary( i )->GetParticle()->GetDefinition()->GetAtomicNumber() == max_SecZ
110  && theResult.Get()->GetSecondary( i )->GetParticle()->GetDefinition()->GetAtomicMass() == max_SecA )
111  {
113  p = p * resi_pd->GetPDGMass()/ G4IonTable::GetIonTable()->GetIon ( max_SecZ , max_SecA , 0.0 )->GetPDGMass();
114  theResult.Get()->GetSecondary( i )->GetParticle()->SetDefinition( resi_pd );
116  }
117  }
118  }
119  }
120  }
121 
122 
123  G4LorentzVector secs_4p_lab( 0.0 );
124 
126  G4double fast = 0;
127  G4double slow = 1;
128  G4int ifast = 0;
129  G4int islow = 0;
130  G4int ires = -1;
131 
132  for ( G4int i = 0 ; i < n_sec ; i++ )
133  {
134 
135  //G4cout << "HP_DB " << i
136  // << " " << theResult.GetSecondary( i )->GetParticle()->GetDefinition()->GetParticleName()
137  // << " 4p " << theResult.GetSecondary( i )->GetParticle()->Get4Momentum()
138  // << " ke " << theResult.GetSecondary( i )->GetParticle()->Get4Momentum().e() - theResult.GetSecondary( i )->GetParticle()->GetDefinition()->GetPDGMass()
139  // << G4endl;
140 
141  secs_4p_lab += theResult.Get()->GetSecondary( i )->GetParticle()->Get4Momentum();
142 
143  G4double beta = 0;
145  {
146  beta = theResult.Get()->GetSecondary( i )->GetParticle()->Get4Momentum().beta();
147  }
148  else
149  {
150  beta = 1;
151  }
152 
153  if ( theResult.Get()->GetSecondary( i )->GetParticle()->GetDefinition() == resi_pd ) ires = i;
154 
155  if ( slow > beta && beta != 0 )
156  {
157  slow = beta;
158  islow = i;
159  }
160 
161  if ( fast <= beta )
162  {
163  if ( fast != 1 )
164  {
165  fast = beta;
166  ifast = i;
167  }
168  else
169  {
170 // fast is already photon then check E
172  if ( e > theResult.Get()->GetSecondary( ifast )->GetParticle()->Get4Momentum().e() )
173  {
174 // among photons, the highest E becomes the fastest
175  ifast = i;
176  }
177  }
178  }
179  }
180 
181 
182  G4LorentzVector dif_4p = init_4p_lab - secs_4p_lab;
183 
184  //G4cout << "HP_DB dif_4p " << init_4p_lab - secs_4p_lab << G4endl;
185  //G4cout << "HP_DB dif_3p mag " << ( dif_4p.v() ).mag() << G4endl;
186  //G4cout << "HP_DB dif_e " << dif_4p.e() - ( dif_4p.v() ).mag()<< G4endl;
187 
188  G4LorentzVector p4(0);
189  if ( ires == -1 )
190  {
191 // Create and Add Residual Nucleus
192  ires = nSecondaries;
193  nSecondaries += 1;
194 
195  G4DynamicParticle* res = new G4DynamicParticle ( resi_pd , dif_4p.v() );
196  theResult.Get()->AddSecondary ( res );
197 
198  p4 = res->Get4Momentum();
199  if ( slow > p4.beta() )
200  {
201  slow = p4.beta();
202  islow = ires;
203  }
204  dif_4p = init_4p_lab - ( secs_4p_lab + p4 );
205  }
206 
207  if ( needOneMoreSec && oneMoreSec_pd)
208  //
209  // fca: this is not a fix, this is a crash avoidance...
210  // fca: the baryon number is still wrong, most probably because it
211  // fca: should have been decreased, but since we could not create a particle
212  // fca: we just do not add it
213  //
214  {
215  nSecondaries += 1;
216  G4DynamicParticle* one = new G4DynamicParticle ( oneMoreSec_pd , dif_4p.v() );
217  theResult.Get()->AddSecondary ( one );
218  p4 = one->Get4Momentum();
219  if ( slow > p4.beta() )
220  {
221  slow = p4.beta();
222  islow = nSecondaries-1; //Because the first is 0th, so the last becomes "nSecondaries-1"
223  }
224  dif_4p = init_4p_lab - ( secs_4p_lab + p4 );
225  }
226 
227  //Which is bigger dif_p or dif_e
228 
229  if ( dif_4p.v().mag() < std::abs( dif_4p.e() ) )
230  {
231 
232  // Adjust p
233  //if ( dif_4p.v().mag() < 1*MeV )
234  if ( minimum_energy < dif_4p.v().mag() && dif_4p.v().mag() < 1*MeV )
235  {
236 
237  nSecondaries += 1;
238  theResult.Get()->AddSecondary ( new G4DynamicParticle ( G4Gamma::Gamma() , dif_4p.v() ) );
239 
240  }
241  else
242  {
243  //G4cout << "HP_DB Difference in dif_p is too large (>1MeV) or too small(<1keV) to adjust, so that give up tuning" << G4endl;
244  }
245 
246  }
247  else
248  {
249 
250  // dif_p > dif_e
251  // at first momentum
252  // Move residual momentum
253 
254  p4 = theResult.Get()->GetSecondary( ires )->GetParticle()->Get4Momentum();
255  theResult.Get()->GetSecondary( ires )->GetParticle()->SetMomentum( p4.v() + dif_4p.v() );
256  dif_4p = init_4p_lab - ( secs_4p_lab - p4 + theResult.Get()->GetSecondary( ires )->GetParticle()->Get4Momentum() );
257 
258  //G4cout << "HP_DB new residual kinetic energy " << theResult.GetSecondary( ires )->GetParticle()->GetKineticEnergy() << G4endl;
259 
260  }
261 
262  G4double dif_e = dif_4p.e() - ( dif_4p.v() ).mag();
263  //G4cout << "HP_DB dif_e " << dif_e << G4endl;
264 
265  if ( dif_e > 0 )
266  {
267 
268 // create 2 gamma
269 
270  nSecondaries += 2;
271  G4double e1 = ( dif_4p.e() -dif_4p.v().mag() ) / 2;
272 
273  if ( minimum_energy < e1 )
274  {
275  G4double costh = 2.*G4UniformRand()-1.;
276  G4double phi = twopi*G4UniformRand();
277  G4ThreeVector dir( std::sin(std::acos(costh))*std::cos(phi),
278  std::sin(std::acos(costh))*std::sin(phi),
279  costh);
280  theResult.Get()->AddSecondary ( new G4DynamicParticle ( G4Gamma::Gamma() , e1*dir ) );
281  theResult.Get()->AddSecondary ( new G4DynamicParticle ( G4Gamma::Gamma() , -e1*dir ) );
282  }
283  else
284  {
285  //G4cout << "HP_DB Difference is too small(<1keV) to adjust, so that neglect it" << G4endl;
286  }
287 
288  }
289  else //dif_e < 0
290  {
291 
292 // At first reduce KE of the fastest secondary;
295  G4ThreeVector dir = ( theResult.Get()->GetSecondary( ifast )->GetParticle()->GetMomentum() ).unit();
296 
297  //G4cout << "HP_DB ifast " << ifast << " ke0 " << ke0 << G4endl;
298 
299  if ( ke0 + dif_e > 0 )
300  {
301  theResult.Get()->GetSecondary( ifast )->GetParticle()->SetKineticEnergy( ke0 + dif_e );
302  G4ThreeVector dp = p0 - theResult.Get()->GetSecondary( ifast )->GetParticle()->GetMomentum();
303 
305  //theResult.GetSecondary( islow )->GetParticle()->SetMomentum( p - dif_e*dir );
306  theResult.Get()->GetSecondary( islow )->GetParticle()->SetMomentum( p + dp );
307  }
308  else
309  {
310  //G4cout << "HP_DB Difference in dif_e too large ( <0MeV ) to adjust, so that give up tuning" << G4endl;
311  }
312 
313  }
314 
315 }
static const double MeV
Definition: G4SIunits.hh:193
void SetMomentum(const G4ThreeVector &momentum)
G4HadSecondary * GetSecondary(size_t i)
G4double GetKineticEnergy() const
G4Cache< G4HadFinalState * > theResult
CLHEP::Hep3Vector G4ThreeVector
static G4NeutronHPManager * GetInstance()
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:463
value_type & Get() const
Definition: G4Cache.hh:282
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
G4int GetAtomicNumber() const
#define G4UniformRand()
Definition: Randomize.hh:93
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
G4int GetAtomicMass() const
void SetKineticEnergy(G4double aEnergy)
G4LorentzVector Get4Momentum() const
static G4IonTable * GetIonTable()
Definition: G4IonTable.hh:78
static const G4double e1
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
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)
void adjust_final_state(G4LorentzVector)
G4int GetNumberOfSecondaries() const
G4ThreeVector GetMomentum() const
CLHEP::HepLorentzVector G4LorentzVector