source: trunk/source/processes/hadronic/models/pre_equilibrium/exciton_model/src/G4PreCompoundFragment.cc @ 1340

Last change on this file since 1340 was 1340, checked in by garnier, 14 years ago

update ti head

File size: 5.1 KB
Line 
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// $Id: G4PreCompoundFragment.cc,v 1.9 2010/08/28 15:16:55 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-03-ref-09 $
28//
29// J. M. Quesada (August 2008). 
30// Based  on previous work by V. Lara
31//
32// Modified:
33// 06.09.2008 JMQ Also external choice has been added for:
34//               - superimposed Coulomb barrier (if useSICB=true)
35// 20.08.2010 V.Ivanchenko cleanup
36//
37
38#include "G4PreCompoundFragment.hh"
39
40G4PreCompoundFragment::
41G4PreCompoundFragment(const G4ParticleDefinition* part,
42                      G4VCoulombBarrier* aCoulombBarrier)
43  : G4VPreCompoundFragment(part,aCoulombBarrier)
44{}
45
46G4PreCompoundFragment::~G4PreCompoundFragment()
47{}
48
49G4double G4PreCompoundFragment::
50CalcEmissionProbability(const G4Fragment & aFragment)
51{
52  //G4cout << theCoulombBarrier << "  " << GetMaximalKineticEnergy() << G4endl;
53  // If  theCoulombBarrier effect is included in the emission probabilities
54  //if (GetMaximalKineticEnergy() <= 0.0)
55  G4double limit = 0.0;
56  if(OPTxs==0 ||  useSICB) { limit = theCoulombBarrier; }
57  if (GetMaximalKineticEnergy() <= limit) 
58    {
59      theEmissionProbability = 0.0;
60      return 0.0;
61    }   
62  // If  theCoulombBarrier effect is included in the emission probabilities
63  //  G4double LowerLimit = 0.;
64  // Coulomb barrier is the lower limit
65  // of integration over kinetic energy
66  G4double LowerLimit = limit;
67
68  // Excitation energy of nucleus after fragment emission is the upper
69  //limit of integration over kinetic energy
70  G4double UpperLimit = GetMaximalKineticEnergy();
71 
72  theEmissionProbability = 
73    IntegrateEmissionProbability(LowerLimit,UpperLimit,aFragment);
74  /*
75  G4cout << "## G4PreCompoundFragment::CalcEmisProb "
76         << "Z= " << aFragment.GetZ_asInt()
77         << " A= " << aFragment.GetA_asInt()
78         << " Elow= " << LowerLimit/MeV
79         << " Eup= " << UpperLimit/MeV
80         << " prob= " << theEmissionProbability
81         << G4endl;
82  */
83  return theEmissionProbability;
84}
85
86G4double G4PreCompoundFragment::
87IntegrateEmissionProbability(G4double Low, G4double Up,
88                             const G4Fragment & aFragment)
89{
90  static const G4int N = 10;
91  // 10-Points Gauss-Legendre abcisas and weights
92  static const G4double w[N] = {
93    0.0666713443086881,
94    0.149451349150581,
95    0.219086362515982,
96    0.269266719309996,
97    0.295524224714753,
98    0.295524224714753,
99    0.269266719309996,
100    0.219086362515982,
101    0.149451349150581,
102    0.0666713443086881
103  };
104  static const G4double x[N] = {
105    -0.973906528517172,
106    -0.865063366688985,
107    -0.679409568299024,
108    -0.433395394129247,
109    -0.148874338981631,
110    0.148874338981631,
111    0.433395394129247,
112    0.679409568299024,
113    0.865063366688985,
114    0.973906528517172
115  };
116 
117  G4double Total = 0.0;
118
119  for (G4int i = 0; i < N; ++i) 
120    {
121      G4double KineticE = 0.5*((Up-Low)*x[i]+(Up+Low));
122      Total += w[i]*ProbabilityDistributionFunction(KineticE, aFragment);
123    }
124  Total  *= 0.5*(Up-Low);
125  return Total;
126}
127
128G4double G4PreCompoundFragment::
129GetKineticEnergy(const G4Fragment & aFragment) 
130{
131  //let's keep this way for consistency with CalcEmissionProbability method
132  G4double V = 0.0;
133  if(OPTxs==0 || useSICB) { V = theCoulombBarrier; }
134
135  G4double Tmax = GetMaximalKineticEnergy(); 
136  if(Tmax < V) { return 0.0; }
137  G4double T(0.0);
138  G4double Probability(1.0);
139  G4double maxProbability = GetEmissionProbability();
140  do 
141    {
142      T = V + G4UniformRand()*(Tmax-V);
143      Probability = ProbabilityDistributionFunction(T,aFragment);     
144    }   while (maxProbability*G4UniformRand() > Probability); 
145  return T;
146}
Note: See TracBrowser for help on using the repository browser.