Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4PionMinusAbsorptionAtRest.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 // G4PionMinusAbsorptionAtRest physics process
27 // Larry Felawka (TRIUMF), April 1998
28 //---------------------------------------------------------------------
29 
30 #include <string.h>
31 #include <cmath>
32 #include <stdio.h>
33 
35 #include "G4DynamicParticle.hh"
36 #include "G4ParticleTypes.hh"
37 #include "G4SystemOfUnits.hh"
38 #include "Randomize.hh"
40  #include "G4HadronicDeprecate.hh"
41 
42 #define MAX_SECONDARIES 100
43 
44 // constructor
45 
46 G4PionMinusAbsorptionAtRest::G4PionMinusAbsorptionAtRest(const G4String& processName,
47  G4ProcessType aType ) :
48  G4VRestProcess (processName, aType), // initialization
49  massPionMinus(G4PionMinus::PionMinus()->GetPDGMass()/GeV),
50  pdefGamma(G4Gamma::Gamma()),
51  pdefPionZero(G4PionZero::PionZero()),
52  pdefPionMinus(G4PionMinus::PionMinus()),
53  pdefProton(G4Proton::Proton()),
54  pdefNeutron(G4Neutron::Neutron()),
55  pdefDeuteron(G4Deuteron::Deuteron()),
56  pdefTriton(G4Triton::Triton()),
57  pdefAlpha(G4Alpha::Alpha())
58 {
59  G4HadronicDeprecate("G4PiMinusAbsorptionAtRest");
60 
61  if (verboseLevel>0) {
62  G4cout << GetProcessName() << " is created "<< G4endl;
63  }
68 
70 }
71 
72 // destructor
73 
75 {
77  delete [] pv;
78  delete [] eve;
79  delete [] gkin;
80 }
81 
83 {
85 }
86 
88 {
90 }
91 
92 // methods.............................................................................
93 
95  const G4ParticleDefinition& particle
96  )
97 {
98  return ( &particle == pdefPionMinus );
99 
100 }
101 
102 // Warning - this method may be optimized away if made "inline"
104 {
105  return ( ngkine );
106 
107 }
108 
109 // Warning - this method may be optimized away if made "inline"
111 {
112  return ( &gkin[0] );
113 
114 }
115 
117  const G4Track& track,
119  )
120 {
121  // beggining of tracking
123 
124  // condition is set to "Not Forced"
125  *condition = NotForced;
126 
127  // get mean life time
128  currentInteractionLength = GetMeanLifeTime(track, condition);
129 
130  if ((currentInteractionLength <0.0) || (verboseLevel>2)){
131  G4cout << "G4PionMinusAbsorptionAtRestProcess::AtRestGetPhysicalInteractionLength ";
132  G4cout << "[ " << GetProcessName() << "]" <<G4endl;
133  track.GetDynamicParticle()->DumpInfo();
134  G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl;
135  G4cout << "MeanLifeTime = " << currentInteractionLength/ns << "[ns]" <<G4endl;
136  }
137 
139 
140 }
141 
143  const G4Track& track,
144  const G4Step&
145  )
146 //
147 // Handles PionMinuss at rest; a PionMinus can either create secondaries or
148 // do nothing (in which case it should be sent back to decay-handling
149 // section
150 //
151 {
152 
153 // Initialize ParticleChange
154 // all members of G4VParticleChange are set to equal to
155 // corresponding member in G4Track
156 
158 
159 // Store some global quantities that depend on current material and particle
160 
161  globalTime = track.GetGlobalTime()/s;
162  G4Material * aMaterial = track.GetMaterial();
163  const G4int numberOfElements = aMaterial->GetNumberOfElements();
164  const G4ElementVector* theElementVector = aMaterial->GetElementVector();
165 
166  const G4double* theAtomicNumberDensity = aMaterial->GetAtomicNumDensityVector();
167  G4double normalization = 0;
168  for ( G4int i1=0; i1 < numberOfElements; i1++ )
169  {
170  normalization += theAtomicNumberDensity[i1] ; // change when nucleon specific
171  // probabilities are included.
172  }
173  G4double runningSum= 0.;
174  G4double random = G4UniformRand()*normalization;
175  for ( G4int i2=0; i2 < numberOfElements; i2++ )
176  {
177  runningSum += theAtomicNumberDensity[i2]; // change when nucleon specific
178  // probabilities are included.
179  if (random<=runningSum)
180  {
181  targetCharge = G4double((*theElementVector)[i2]->GetZ());
182  targetAtomicMass = (*theElementVector)[i2]->GetN();
183  }
184  }
185  if (random>runningSum)
186  {
187  targetCharge = G4double((*theElementVector)[numberOfElements-1]->GetZ());
188  targetAtomicMass = (*theElementVector)[numberOfElements-1]->GetN();
189 
190  }
191 
192  if (verboseLevel>1) {
193  G4cout << "G4PionMinusAbsorptionAtRest::AtRestDoIt is invoked " <<G4endl;
194  }
195 
196  G4ParticleMomentum momentum;
197  G4float localtime;
198 
200 
201  GenerateSecondaries(); // Generate secondaries
202 
204 
205  for ( G4int isec = 0; isec < ngkine; isec++ ) {
206  G4DynamicParticle* aNewParticle = new G4DynamicParticle;
207  aNewParticle->SetDefinition( gkin[isec].GetParticleDef() );
208  aNewParticle->SetMomentum( gkin[isec].GetMomentum() * GeV );
209 
210  localtime = globalTime + gkin[isec].GetTOF();
211 
212  G4Track* aNewTrack = new G4Track( aNewParticle, localtime*s, position );
213  aNewTrack->SetTouchableHandle(track.GetTouchableHandle());
214  aParticleChange.AddSecondary( aNewTrack );
215 
216  }
217 
219 
220  aParticleChange.ProposeTrackStatus(fStopAndKill); // Kill the incident PionMinus
221 
222 // clear InteractionLengthLeft
223 
225 
226  return &aParticleChange;
227 
228 }
229 
230 
231 void G4PionMinusAbsorptionAtRest::GenerateSecondaries()
232 {
233  static G4int index;
234  static G4int l;
235  static G4int nopt;
236  static G4int i;
237  // DHW 15 May 2011: unused: static G4ParticleDefinition* jnd;
238 
239  for (i = 1; i <= MAX_SECONDARIES; ++i) {
240  pv[i].SetZero();
241  }
242 
243  ngkine = 0; // number of generated secondary particles
244  ntot = 0;
245  result.SetZero();
246  result.SetMass( massPionMinus );
247  result.SetKineticEnergyAndUpdate( 0. );
248  result.SetTOF( 0. );
249  result.SetParticleDef( pdefPionMinus );
250 
251  PionMinusAbsorption(&nopt);
252 
253  // *** CHECK WHETHER THERE ARE NEW PARTICLES GENERATED ***
254  if (ntot != 0 || result.GetParticleDef() != pdefPionMinus) {
255  // *** CURRENT PARTICLE IS NOT THE SAME AS IN THE BEGINNING OR/AND ***
256  // *** ONE OR MORE SECONDARIES HAVE BEEN GENERATED ***
257 
258  // --- INITIAL PARTICLE TYPE HAS BEEN CHANGED ==> PUT NEW TYPE ON ---
259  // --- THE GEANT TEMPORARY STACK ---
260 
261  // --- PUT PARTICLE ON THE STACK ---
262  gkin[0] = result;
263  gkin[0].SetTOF( result.GetTOF() * 5e-11 );
264  ngkine = 1;
265 
266  // --- ALL QUANTITIES ARE TAKEN FROM THE GHEISHA STACK WHERE THE ---
267  // --- CONVENTION IS THE FOLLOWING ---
268 
269  // --- ONE OR MORE SECONDARIES HAVE BEEN GENERATED ---
270  for (l = 1; l <= ntot; ++l) {
271  index = l - 1;
272  // DHW 15 May 2011: unused: jnd = eve[index].GetParticleDef();
273 
274  // --- ADD PARTICLE TO THE STACK IF STACK NOT YET FULL ---
275  if (ngkine < MAX_SECONDARIES) {
276  gkin[ngkine] = eve[index];
277  gkin[ngkine].SetTOF( eve[index].GetTOF() * 5e-11 );
278  ++ngkine;
279  }
280  }
281  }
282  else {
283  // --- NO SECONDARIES GENERATED AND PARTICLE IS STILL THE SAME ---
284  // --- ==> COPY EVERYTHING BACK IN THE CURRENT GEANT STACK ---
285  ngkine = 0;
286  ntot = 0;
287  globalTime += result.GetTOF() * G4float(5e-11);
288  }
289 
290  // --- LIMIT THE VALUE OF NGKINE IN CASE OF OVERFLOW ---
291  ngkine = G4int(std::min(ngkine,G4int(MAX_SECONDARIES)));
292 
293 } // GenerateSecondaries
294 
295 
296 void G4PionMinusAbsorptionAtRest::PionMinusAbsorption(G4int *nopt)
297 {
298  static G4int i;
299  static G4int nt, nbl;
300  static G4float ran, tex;
301  static G4int isw;
302  static G4float ran2, tof1, ekin;
303  static G4float ekin1, ekin2, black;
304  static G4float pnrat;
305  static G4ParticleDefinition* ipa1;
306  static G4ParticleDefinition* inve;
307 
308  // *** CHARGED PION ABSORPTION BY A NUCLEUS ***
309  // *** NVE 04-MAR-1988 CERN GENEVA ***
310 
311  // ORIGIN : H.FESEFELDT (09-JULY-1987)
312 
313  // PANOFSKY RATIO (PI- P --> N PI0/PI- P --> N GAMMA) = 3/2
314  // FOR CAPTURE ON PROTON (HYDROGEN),
315  // STAR PRODUCTION FOR HEAVIER ELEMENTS
316 
317  pv[1].SetZero();
318  pv[1].SetMass( massPionMinus );
319  pv[1].SetKineticEnergyAndUpdate( 0. );
320  pv[1].SetTOF( result.GetTOF() );
321  pv[1].SetParticleDef( result.GetParticleDef() );
322  if (targetAtomicMass <= G4float(1.5)) {
323  ran = G4UniformRand();
324  isw = 1;
325  if (ran < G4float(.33)) {
326  isw = 2;
327  }
328  *nopt = isw;
329  ran = G4UniformRand();
330  tof1 = std::log(ran) * G4float(-25.);
331  tof1 *= G4float(20.);
332  if (isw != 1) {
333  pv[2].SetZero();
334  pv[2].SetMass( 0. );
335  pv[2].SetKineticEnergyAndUpdate( .02 );
336  pv[2].SetTOF( result.GetTOF() + tof1 );
337  pv[2].SetParticleDef( pdefGamma );
338  }
339  else {
340  pv[2] = pv[1];
341  pv[2].SetTOF( result.GetTOF() + tof1 );
342  pv[2].SetParticleDef( pdefPionZero );
343  }
344  result = pv[2];
345  }
346  else {
347  // **
348  // ** STAR PRODUCTION FOR PION ABSORPTION IN HEAVY ELEMENTS
349  // **
350  evapEnergy1 = G4float(.0135);
351  evapEnergy3 = G4float(.0058);
352  nt = 1;
353  tex = evapEnergy1;
354  black = std::log(targetAtomicMass) * G4float(.5);
355  Poisso(black, &nbl);
356  if (nbl <= 0) {
357  nbl = 1;
358  }
359  if (nt + nbl > (MAX_SECONDARIES - 2)) {
360  nbl = (MAX_SECONDARIES - 2) - nt;
361  }
362  ekin = tex / nbl;
363  ekin2 = G4float(0.);
364  for (i = 1; i <= nbl; ++i) {
365  if (nt == (MAX_SECONDARIES - 2)) {
366  continue;
367  }
368  ran2 = G4UniformRand();
369  ekin1 = -G4double(ekin) * std::log(ran2);
370  ekin2 += ekin1;
371  ipa1 = pdefNeutron;
372  pnrat = G4float(1.) - targetCharge / targetAtomicMass;
373  if (G4UniformRand() > pnrat) {
374  ipa1 = pdefProton;
375  }
376  ++nt;
377  pv[nt].SetZero();
378  pv[nt].SetMass( ipa1->GetPDGMass()/GeV );
379  pv[nt].SetKineticEnergyAndUpdate( ekin1 );
380  pv[nt].SetTOF( 2. );
381  pv[nt].SetParticleDef( ipa1 );
382  if (ekin2 > tex) {
383  break;
384  }
385  }
386  tex = evapEnergy3;
387  black = std::log(targetAtomicMass) * G4float(.5);
388  Poisso(black, &nbl);
389  if (nt + nbl > (MAX_SECONDARIES - 2)) {
390  nbl = (MAX_SECONDARIES - 2) - nt;
391  }
392  if (nbl <= 0) {
393  nbl = 1;
394  }
395  ekin = tex / nbl;
396  ekin2 = G4float(0.);
397  for (i = 1; i <= nbl; ++i) {
398  if (nt == (MAX_SECONDARIES - 2)) {
399  continue;
400  }
401  ran2 = G4UniformRand();
402  ekin1 = -G4double(ekin) * std::log(ran2);
403  ekin2 += ekin1;
404  ++nt;
405  ran = G4UniformRand();
406  inve= pdefDeuteron;
407  if (ran > G4float(.6)) {
408  inve = pdefTriton;
409  }
410  if (ran > G4float(.9)) {
411  inve = pdefAlpha;
412  }
413  pv[nt].SetZero();
414  pv[nt].SetMass( inve->GetPDGMass()/GeV );
415  pv[nt].SetKineticEnergyAndUpdate( ekin1 );
416  pv[nt].SetTOF( 2. );
417  pv[nt].SetParticleDef( inve );
418  if (ekin2 > tex) {
419  break;
420  }
421  }
422  // **
423  // ** STORE ON EVENT COMMON
424  // **
425  ran = G4UniformRand();
426  tof1 = std::log(ran) * G4float(-25.);
427  tof1 *= G4float(20.);
428  for (i = 2; i <= nt; ++i) {
429  pv[i].SetTOF( result.GetTOF() + tof1 );
430  }
431  result = pv[2];
432  for (i = 3; i <= nt; ++i) {
433  if (ntot >= MAX_SECONDARIES) {
434  break;
435  }
436  eve[ntot++] = pv[i];
437  }
438  }
439 
440 } // PionMinusAbsorption
441 
442 
443 void G4PionMinusAbsorptionAtRest::Poisso(G4float xav, G4int *iran)
444 {
445  static G4int i;
446  static G4float r, p1, p2, p3;
447  static G4int fivex;
448  static G4float rr, ran, rrr, ran1;
449 
450  // *** GENERATION OF POISSON DISTRIBUTION ***
451  // *** NVE 16-MAR-1988 CERN GENEVA ***
452  // ORIGIN : H.FESEFELDT (27-OCT-1983)
453 
454  // --- USE NORMAL DISTRIBUTION FOR <X> > 9.9 ---
455  if (xav > G4float(9.9)) {
456  // ** NORMAL DISTRIBUTION WITH SIGMA**2 = <X>
457  Normal(&ran1);
458  ran1 = xav + ran1 * std::sqrt(xav);
459  *iran = G4int(ran1);
460  if (*iran < 0) {
461  *iran = 0;
462  }
463  }
464  else {
465  fivex = G4int(xav * G4float(5.));
466  *iran = 0;
467  if (fivex > 0) {
468  r = std::exp(-G4double(xav));
469  ran1 = G4UniformRand();
470  if (ran1 > r) {
471  rr = r;
472  for (i = 1; i <= fivex; ++i) {
473  ++(*iran);
474  if (i <= 5) {
475  rrr = std::pow(xav, G4float(i)) / NFac(i);
476  }
477  // ** STIRLING' S FORMULA FOR LARGE NUMBERS
478  if (i > 5) {
479  rrr = std::exp(i * std::log(xav) -
480  (i + G4float(.5)) * std::log(i * G4float(1.)) +
481  i - G4float(.9189385));
482  }
483  rr += r * rrr;
484  if (ran1 <= rr) {
485  break;
486  }
487  }
488  }
489  }
490  else {
491  // ** FOR VERY SMALL XAV TRY IRAN=1,2,3
492  p1 = xav * std::exp(-G4double(xav));
493  p2 = xav * p1 / G4float(2.);
494  p3 = xav * p2 / G4float(3.);
495  ran = G4UniformRand();
496  if (ran >= p3) {
497  if (ran >= p2) {
498  if (ran >= p1) {
499  *iran = 0;
500  }
501  else {
502  *iran = 1;
503  }
504  }
505  else {
506  *iran = 2;
507  }
508  }
509  else {
510  *iran = 3;
511  }
512  }
513  }
514 
515 } // Poisso
516 
517 
518 G4int G4PionMinusAbsorptionAtRest::NFac(G4int n)
519 {
520  G4int ret_val;
521  static G4int i, j;
522 
523  // *** NVE 16-MAR-1988 CERN GENEVA ***
524  // ORIGIN : H.FESEFELDT (27-OCT-1983)
525 
526  ret_val = 1;
527  j = n;
528  if (j > 1) {
529  if (j > 10) {
530  j = 10;
531  }
532  for (i = 2; i <= j; ++i) {
533  ret_val *= i;
534  }
535  }
536  return ret_val;
537 
538 } // NFac
539 
540 
541 void G4PionMinusAbsorptionAtRest::Normal(G4float *ran)
542 {
543  static G4int i;
544 
545  // *** NVE 14-APR-1988 CERN GENEVA ***
546  // ORIGIN : H.FESEFELDT (27-OCT-1983)
547 
548  *ran = G4float(-6.);
549  for (i = 1; i <= 12; ++i) {
550  *ran += G4UniformRand();
551  }
552 
553 } // Normal