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

Protected Member Functions

 G4EvaporationChannel ()

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 51 of file G4EvaporationChannel.cc.

References G4NucleiProperties::GetNuclearMass().

00054                                                                               :
00055     G4VEvaporationChannel(aName),
00056     theA(anA),
00057     theZ(aZ),
00058     theEvaporationProbabilityPtr(aEmissionStrategy),
00059     theCoulombBarrierPtr(aCoulombBarrier),
00060     EmissionProbability(0.0),
00061     MaximalKineticEnergy(-1000.0)
00062 { 
00063   ResidualA = 0;
00064   ResidualZ = 0;
00065   ResidualMass = CoulombBarrier=0.0;
00066   EvaporatedMass = G4NucleiProperties::GetNuclearMass(theA, theZ);
00067   theLevelDensityPtr = new G4EvaporationLevelDensityParameter;
00068 }

G4EvaporationChannel::~G4EvaporationChannel (  )  [virtual]

Definition at line 85 of file G4EvaporationChannel.cc.

00086 {
00087   delete theLevelDensityPtr;
00088 }

G4EvaporationChannel::G4EvaporationChannel (  )  [protected]

Definition at line 70 of file G4EvaporationChannel.cc.

00070                                           :
00071     G4VEvaporationChannel(""),
00072     theA(0),
00073     theZ(0),
00074     theEvaporationProbabilityPtr(0),
00075     theCoulombBarrierPtr(0),
00076     EmissionProbability(0.0),
00077     MaximalKineticEnergy(-1000.0)
00078 { 
00079   ResidualA = 0;
00080   ResidualZ = 0;
00081   EvaporatedMass = ResidualMass = CoulombBarrier = 0.0;
00082   theLevelDensityPtr = new G4EvaporationLevelDensityParameter;
00083 }


Member Function Documentation

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

Implements G4VEvaporationChannel.

Definition at line 151 of file G4EvaporationChannel.cc.

References G4cout, G4endl, and G4Fragment::GetMomentum().

00152 {
00153   /*
00154   G4double Ecm = GetKineticEnergy(theNucleus) + ResidualMass + EvaporatedMass;
00155   
00156   G4double EvaporatedEnergy = 
00157     ((Ecm-ResidualMass)*(Ecm+ResidualMass) + EvaporatedMass*EvaporatedMass)/(2*Ecm);
00158   */
00159   G4double EvaporatedEnergy = GetKineticEnergy(theNucleus) + EvaporatedMass;
00160 
00161   G4ThreeVector momentum(IsotropicVector
00162                          (std::sqrt((EvaporatedEnergy - EvaporatedMass)*
00163                                     (EvaporatedEnergy + EvaporatedMass))));
00164   
00165   G4LorentzVector EvaporatedMomentum(momentum,EvaporatedEnergy);
00166   G4LorentzVector ResidualMomentum = theNucleus.GetMomentum();
00167   EvaporatedMomentum.boost(ResidualMomentum.boostVector());
00168   
00169   G4Fragment * EvaporatedFragment = new G4Fragment(theA,theZ,EvaporatedMomentum);
00170   ResidualMomentum -= EvaporatedMomentum;
00171 
00172   G4Fragment * ResidualFragment = new G4Fragment(ResidualA, ResidualZ, ResidualMomentum);
00173  
00174   G4FragmentVector * theResult = new G4FragmentVector;
00175   
00176 #ifdef debug
00177   G4double Efinal = ResidualMomentum.e() + EvaporatedMomentum.e();
00178   G4ThreeVector Pfinal = ResidualMomentum.vect() + EvaporatedMomentum.vect();
00179   if (std::abs(Efinal-theNucleus.GetMomentum().e()) > 1.0*keV) {
00180     G4cout << "@@@@@@@@@@@@@@@@@@@@@ G4Evaporation Chanel: ENERGY @@@@@@@@@@@@@@@@" << G4endl;
00181     G4cout << "Initial : " << theNucleus.GetMomentum().e()/MeV << " MeV    Final :" 
00182            <<Efinal/MeV << " MeV    Delta: " <<  (Efinal-theNucleus.GetMomentum().e())/MeV 
00183            << " MeV" << G4endl;
00184   }
00185   if (std::abs(Pfinal.x()-theNucleus.GetMomentum().x()) > 1.0*keV ||
00186       std::abs(Pfinal.y()-theNucleus.GetMomentum().y()) > 1.0*keV ||
00187       std::abs(Pfinal.z()-theNucleus.GetMomentum().z()) > 1.0*keV ) {
00188     G4cout << "@@@@@@@@@@@@@@@@@@@@@ G4Evaporation Chanel: MOMENTUM @@@@@@@@@@@@@@@@" << G4endl;
00189     G4cout << "Initial : " << theNucleus.GetMomentum().vect() << " MeV    Final :" 
00190            <<Pfinal/MeV << " MeV    Delta: " <<  Pfinal-theNucleus.GetMomentum().vect()
00191            << " MeV" << G4endl;   
00192     
00193   }
00194 #endif
00195   theResult->push_back(EvaporatedFragment);
00196   theResult->push_back(ResidualFragment);
00197   return theResult; 
00198 } 

G4double G4EvaporationChannel::GetEmissionProbability ( G4Fragment fragment  )  [virtual]

Implements G4VEvaporationChannel.

Definition at line 93 of file G4EvaporationChannel.cc.

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

00094 {
00095   //for inverse cross section choice
00096   theEvaporationProbabilityPtr->SetOPTxs(OPTxs);
00097   // for superimposed Coulomb Barrier for inverse cross sections
00098   theEvaporationProbabilityPtr->UseSICB(useSICB);
00099   
00100   G4int FragmentA = fragment->GetA_asInt();
00101   G4int FragmentZ = fragment->GetZ_asInt();
00102   ResidualA = FragmentA - theA;
00103   ResidualZ = FragmentZ - theZ;
00104   //G4cout << "G4EvaporationChannel::Initialize Z= " << theZ << " A= " << theA 
00105   //     << " FragZ= " << FragmentZ << " FragA= " << FragmentA << G4endl;
00106   
00107   //Effective excitation energy
00108   G4double ExEnergy = fragment->GetExcitationEnergy() - 
00109     G4PairingCorrection::GetInstance()->GetPairingCorrection(FragmentA,FragmentZ);
00110   
00111   // Only channels which are physically allowed are taken into account 
00112   if (ResidualA <= 0 || ResidualZ <= 0 || ResidualA < ResidualZ ||
00113       (ResidualA == ResidualZ && ResidualA > 1) || ExEnergy <= 0.0) {
00114     CoulombBarrier = ResidualMass = 0.0;
00115     MaximalKineticEnergy = -1000.0*MeV;
00116     EmissionProbability = 0.0;
00117   } else {
00118     ResidualMass = G4NucleiProperties::GetNuclearMass(ResidualA, ResidualZ);
00119     G4double FragmentMass = fragment->GetGroundStateMass();
00120     CoulombBarrier = theCoulombBarrierPtr->GetCoulombBarrier(ResidualA,ResidualZ,ExEnergy);
00121     // Maximal Kinetic Energy
00122     MaximalKineticEnergy = CalcMaximalKineticEnergy(FragmentMass + ExEnergy);
00123     //MaximalKineticEnergy = ExEnergy + fragment.GetGroundStateMass() 
00124     //  - EvaporatedMass - ResidualMass;
00125     
00126     // Emission probability
00127     // Protection for the case Tmax<V. If not set in this way we could end up in an 
00128     // infinite loop in  the method GetKineticEnergy if OPTxs!=0 && useSICB=true. 
00129     // Of course for OPTxs=0 we have the Coulomb barrier 
00130     
00131     G4double limit;
00132     if (OPTxs==0 || (OPTxs!=0 && useSICB)) 
00133       limit= CoulombBarrier;
00134     else limit=0.;
00135     limit = std::max(limit, FragmentMass - ResidualMass - EvaporatedMass);
00136   
00137     // The threshold for charged particle emission must be  set to 0 if Coulomb 
00138     //cutoff  is included in the cross sections
00139     if (MaximalKineticEnergy <= limit) EmissionProbability = 0.0;
00140     else { 
00141       // Total emission probability for this channel
00142       EmissionProbability = theEvaporationProbabilityPtr->
00143         EmissionProbability(*fragment, MaximalKineticEnergy);
00144     }
00145   }
00146   //G4cout << "G4EvaporationChannel:: probability= " << EmissionProbability << G4endl; 
00147   
00148   return EmissionProbability;
00149 }

G4double G4EvaporationChannel::GetMaximalKineticEnergy ( void   )  const [inline]

Definition at line 76 of file G4EvaporationChannel.hh.

00077   { return MaximalKineticEnergy; }

void G4EvaporationChannel::SetCoulombBarrierStrategy ( G4VCoulombBarrier aCoulombBarrier  )  [inline]

Definition at line 61 of file G4EvaporationChannel.hh.

00062   {theCoulombBarrierPtr = aCoulombBarrier;} 

void G4EvaporationChannel::SetEmissionStrategy ( G4EvaporationProbability aEmissionStrategy  )  [inline]

Definition at line 58 of file G4EvaporationChannel.hh.

00059   {theEvaporationProbabilityPtr = aEmissionStrategy;}


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:51:56 2013 for Geant4 by  doxygen 1.4.7