Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Protected Member Functions
G4EvaporationChannel Class Reference

#include <G4EvaporationChannel.hh>

Inheritance diagram for G4EvaporationChannel:
G4VEvaporationChannel G4AlphaEvaporationChannel G4DeuteronEvaporationChannel G4He3EvaporationChannel G4NeutronEvaporationChannel G4ProtonEvaporationChannel G4TritonEvaporationChannel

Public Member Functions

 G4EvaporationChannel (G4int theA, G4int theZ, const G4String &aName, G4EvaporationProbability *aEmissionStrategy, G4VCoulombBarrier *aCoulombBarrier)
 
virtual ~G4EvaporationChannel ()
 
void SetEmissionStrategy (G4EvaporationProbability *aEmissionStrategy)
 
void SetCoulombBarrierStrategy (G4VCoulombBarrier *aCoulombBarrier)
 
virtual G4double GetEmissionProbability (G4Fragment *fragment)
 
G4FragmentVectorBreakUp (const G4Fragment &theNucleus)
 
G4double GetMaximalKineticEnergy (void) const
 
- Public Member Functions inherited from G4VEvaporationChannel
 G4VEvaporationChannel (const G4String &aName="Anonymous", G4EvaporationChannelType timeType=fDelayedEmission)
 
virtual ~G4VEvaporationChannel ()
 
virtual G4double GetLifeTime (G4Fragment *theNucleus)
 
virtual G4FragmentEmittedFragment (G4Fragment *theNucleus)
 
virtual G4FragmentVectorBreakUpFragment (G4Fragment *theNucleus)
 
G4String GetName () const
 
void SetName (const G4String &aName)
 
void SetOPTxs (G4int opt)
 
void UseSICB (G4bool use)
 

Protected Member Functions

 G4EvaporationChannel ()
 

Additional Inherited Members

- Protected Attributes inherited from G4VEvaporationChannel
G4EvaporationChannelType sampleDecayTime
 
G4int OPTxs
 
G4bool useSICB
 

Detailed Description

Definition at line 47 of file G4EvaporationChannel.hh.

Constructor & Destructor Documentation

G4EvaporationChannel::G4EvaporationChannel ( G4int  theA,
G4int  theZ,
const G4String aName,
G4EvaporationProbability aEmissionStrategy,
G4VCoulombBarrier aCoulombBarrier 
)

Definition at line 53 of file G4EvaporationChannel.cc.

References G4NucleiProperties::GetNuclearMass().

56  :
57  G4VEvaporationChannel(aName),
58  theA(anA),
59  theZ(aZ),
60  theEvaporationProbabilityPtr(aEmissionStrategy),
61  theCoulombBarrierPtr(aCoulombBarrier),
62  EmissionProbability(0.0),
63  MaximalKineticEnergy(-1000.0)
64 {
65  ResidualA = 0;
66  ResidualZ = 0;
67  ResidualMass = CoulombBarrier=0.0;
68  EvaporatedMass = G4NucleiProperties::GetNuclearMass(theA, theZ);
69  theLevelDensityPtr = new G4EvaporationLevelDensityParameter;
70 }
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4VEvaporationChannel(const G4String &aName="Anonymous", G4EvaporationChannelType timeType=fDelayedEmission)
G4EvaporationChannel::~G4EvaporationChannel ( )
virtual

Definition at line 87 of file G4EvaporationChannel.cc.

88 {
89  delete theLevelDensityPtr;
90 }
G4EvaporationChannel::G4EvaporationChannel ( )
protected

Definition at line 72 of file G4EvaporationChannel.cc.

72  :
74  theA(0),
75  theZ(0),
76  theEvaporationProbabilityPtr(0),
77  theCoulombBarrierPtr(0),
78  EmissionProbability(0.0),
79  MaximalKineticEnergy(-1000.0)
80 {
81  ResidualA = 0;
82  ResidualZ = 0;
83  EvaporatedMass = ResidualMass = CoulombBarrier = 0.0;
84  theLevelDensityPtr = new G4EvaporationLevelDensityParameter;
85 }
G4VEvaporationChannel(const G4String &aName="Anonymous", G4EvaporationChannelType timeType=fDelayedEmission)

Member Function Documentation

G4FragmentVector * G4EvaporationChannel::BreakUp ( const G4Fragment theNucleus)
virtual

Implements G4VEvaporationChannel.

Definition at line 154 of file G4EvaporationChannel.cc.

References CLHEP::HepLorentzVector::boost(), CLHEP::HepLorentzVector::boostVector(), CLHEP::HepLorentzVector::e(), G4cout, G4endl, G4Fragment::GetMomentum(), python.hepunit::keV, python.hepunit::MeV, CLHEP::HepLorentzVector::vect(), CLHEP::Hep3Vector::x(), CLHEP::HepLorentzVector::x(), CLHEP::Hep3Vector::y(), CLHEP::HepLorentzVector::y(), CLHEP::Hep3Vector::z(), and CLHEP::HepLorentzVector::z().

155 {
156  /*
157  G4double Ecm = GetKineticEnergy(theNucleus) + ResidualMass + EvaporatedMass;
158 
159  G4double EvaporatedEnergy =
160  ((Ecm-ResidualMass)*(Ecm+ResidualMass) + EvaporatedMass*EvaporatedMass)/(2*Ecm);
161  */
162  G4double EvaporatedEnergy = GetKineticEnergy(theNucleus) + EvaporatedMass;
163 
164  G4ThreeVector momentum(IsotropicVector
165  (std::sqrt((EvaporatedEnergy - EvaporatedMass)*
166  (EvaporatedEnergy + EvaporatedMass))));
167 
168  G4LorentzVector EvaporatedMomentum(momentum,EvaporatedEnergy);
169  G4LorentzVector ResidualMomentum = theNucleus.GetMomentum();
170  EvaporatedMomentum.boost(ResidualMomentum.boostVector());
171 
172  G4Fragment * EvaporatedFragment = new G4Fragment(theA,theZ,EvaporatedMomentum);
173  ResidualMomentum -= EvaporatedMomentum;
174 
175  G4Fragment * ResidualFragment = new G4Fragment(ResidualA, ResidualZ, ResidualMomentum);
176 
177  G4FragmentVector * theResult = new G4FragmentVector;
178 
179 #ifdef debug
180  G4double Efinal = ResidualMomentum.e() + EvaporatedMomentum.e();
181  G4ThreeVector Pfinal = ResidualMomentum.vect() + EvaporatedMomentum.vect();
182  if (std::abs(Efinal-theNucleus.GetMomentum().e()) > 1.0*keV) {
183  G4cout << "@@@@@@@@@@@@@@@@@@@@@ G4Evaporation Chanel: ENERGY @@@@@@@@@@@@@@@@"
184  << G4endl;
185  G4cout << "Initial : " << theNucleus.GetMomentum().e()/MeV << " MeV Final :"
186  <<Efinal/MeV << " MeV Delta: "
187  << (Efinal-theNucleus.GetMomentum().e())/MeV
188  << " MeV" << G4endl;
189  }
190  if (std::abs(Pfinal.x()-theNucleus.GetMomentum().x()) > 1.0*keV ||
191  std::abs(Pfinal.y()-theNucleus.GetMomentum().y()) > 1.0*keV ||
192  std::abs(Pfinal.z()-theNucleus.GetMomentum().z()) > 1.0*keV ) {
193  G4cout << "@@@@@@@@@@@@@@@@@@@@@ G4Evaporation Chanel: MOMENTUM @@@@@@@@@@@@@@@@"
194  << G4endl;
195  G4cout << "Initial : " << theNucleus.GetMomentum().vect() << " MeV Final :"
196  <<Pfinal/MeV << " MeV Delta: " << Pfinal-theNucleus.GetMomentum().vect()
197  << " MeV" << G4endl;
198 
199  }
200 #endif
201  theResult->push_back(EvaporatedFragment);
202  theResult->push_back(ResidualFragment);
203  return theResult;
204 }
Hep3Vector boostVector() const
double x() const
double z() const
Hep3Vector vect() const
G4GLOB_DLL std::ostream G4cout
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:271
HepLorentzVector & boost(double, double, double)
std::vector< G4Fragment * > G4FragmentVector
Definition: G4Fragment.hh:65
double y() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4double G4EvaporationChannel::GetEmissionProbability ( G4Fragment fragment)
virtual

Implements G4VEvaporationChannel.

Definition at line 95 of file G4EvaporationChannel.cc.

References G4Fragment::GetA_asInt(), G4VCoulombBarrier::GetCoulombBarrier(), G4Fragment::GetExcitationEnergy(), G4Fragment::GetGroundStateMass(), G4PairingCorrection::GetInstance(), G4NucleiProperties::GetNuclearMass(), G4PairingCorrection::GetPairingCorrection(), G4Fragment::GetZ_asInt(), G4INCL::Math::max(), python.hepunit::MeV, G4VEvaporationChannel::OPTxs, G4VEmissionProbability::SetOPTxs(), G4VEmissionProbability::UseSICB(), and G4VEvaporationChannel::useSICB.

96 {
97  //for inverse cross section choice
98  theEvaporationProbabilityPtr->SetOPTxs(OPTxs);
99  // for superimposed Coulomb Barrier for inverse cross sections
100  theEvaporationProbabilityPtr->UseSICB(useSICB);
101 
102  G4int FragmentA = fragment->GetA_asInt();
103  G4int FragmentZ = fragment->GetZ_asInt();
104  ResidualA = FragmentA - theA;
105  ResidualZ = FragmentZ - theZ;
106  //G4cout << "G4EvaporationChannel::Initialize Z= " << theZ << " A= " << theA
107  // << " FragZ= " << FragmentZ << " FragA= " << FragmentA << G4endl;
108 
109  //Effective excitation energy
110  G4double ExEnergy = fragment->GetExcitationEnergy() -
112 
113  // Only channels which are physically allowed are taken into account
114  if (ResidualA <= 0 || ResidualZ <= 0 || ResidualA < ResidualZ ||
115  (ResidualA == ResidualZ && ResidualA > 1) || ExEnergy <= 0.0) {
116  CoulombBarrier = ResidualMass = 0.0;
117  MaximalKineticEnergy = -1000.0*MeV;
118  EmissionProbability = 0.0;
119  } else {
120  ResidualMass = G4NucleiProperties::GetNuclearMass(ResidualA, ResidualZ);
121  G4double FragmentMass = fragment->GetGroundStateMass();
122  CoulombBarrier =
123  theCoulombBarrierPtr->GetCoulombBarrier(ResidualA,ResidualZ,ExEnergy);
124  // Maximal Kinetic Energy
125  MaximalKineticEnergy = CalcMaximalKineticEnergy(FragmentMass + ExEnergy);
126  //MaximalKineticEnergy = ExEnergy + fragment.GetGroundStateMass()
127  // - EvaporatedMass - ResidualMass;
128 
129  // Emission probability
130  // Protection for the case Tmax<V. If not set in this way we could end up in an
131  // infinite loop in the method GetKineticEnergy if OPTxs!=0 && useSICB=true.
132  // Of course for OPTxs=0 we have the Coulomb barrier
133 
134  G4double limit;
135  if (OPTxs==0 || (OPTxs!=0 && useSICB))
136  limit= CoulombBarrier;
137  else limit=0.;
138  limit = std::max(limit, FragmentMass - ResidualMass - EvaporatedMass);
139 
140  // The threshold for charged particle emission must be set to 0 if Coulomb
141  //cutoff is included in the cross sections
142  if (MaximalKineticEnergy <= limit) EmissionProbability = 0.0;
143  else {
144  // Total emission probability for this channel
145  EmissionProbability = theEvaporationProbabilityPtr->
146  EmissionProbability(*fragment, MaximalKineticEnergy);
147  }
148  }
149  //G4cout << "G4EvaporationChannel:: probability= " << EmissionProbability << G4endl;
150 
151  return EmissionProbability;
152 }
static G4double GetNuclearMass(const G4double A, const G4double Z)
int G4int
Definition: G4Types.hh:78
G4int GetA_asInt() const
Definition: G4Fragment.hh:238
G4double GetPairingCorrection(G4int A, G4int Z) const
G4double GetGroundStateMass() const
Definition: G4Fragment.hh:260
static G4PairingCorrection * GetInstance()
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4int GetZ_asInt() const
Definition: G4Fragment.hh:243
virtual G4double GetCoulombBarrier(G4int ARes, G4int ZRes, G4double U) const =0
double G4double
Definition: G4Types.hh:76
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:255
G4double G4EvaporationChannel::GetMaximalKineticEnergy ( void  ) const
inline

Definition at line 76 of file G4EvaporationChannel.hh.

77  { return MaximalKineticEnergy; }
void G4EvaporationChannel::SetCoulombBarrierStrategy ( G4VCoulombBarrier aCoulombBarrier)
inline

Definition at line 61 of file G4EvaporationChannel.hh.

62  {theCoulombBarrierPtr = aCoulombBarrier;}
void G4EvaporationChannel::SetEmissionStrategy ( G4EvaporationProbability aEmissionStrategy)
inline

Definition at line 58 of file G4EvaporationChannel.hh.

59  {theEvaporationProbabilityPtr = aEmissionStrategy;}

The documentation for this class was generated from the following files: