G4PreCompoundEmission.cc

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 // $Id$
00027 //
00028 // -------------------------------------------------------------------
00029 //
00030 // GEANT4 Class file
00031 //
00032 //
00033 // File name:     G4PreCompoundEmission
00034 //
00035 // Author:         V.Lara
00036 //
00037 // Modified:  
00038 // 15.01.2010 J.M.Quesada  added protection against unphysical values of parameter an 
00039 // 19.01.2010 V.Ivanchenko simplified computation of parameter an, sample cosTheta 
00040 //                         instead of theta; protect all calls to sqrt 
00041 // 20.08.2010 V.Ivanchenko added G4Pow and G4PreCompoundParameters pointers
00042 //                         use int Z and A and cleanup
00043 //
00044 
00045 #include "G4PreCompoundEmission.hh"
00046 #include "G4PhysicalConstants.hh"
00047 #include "G4SystemOfUnits.hh"
00048 #include "G4Pow.hh"
00049 #include "Randomize.hh"
00050 #include "G4PreCompoundParameters.hh"
00051 #include "G4PreCompoundEmissionFactory.hh"
00052 #include "G4HETCEmissionFactory.hh"
00053 #include "G4HadronicException.hh"
00054 
00055 G4PreCompoundEmission::G4PreCompoundEmission()
00056 {
00057   theFragmentsFactory = new G4PreCompoundEmissionFactory();
00058   theFragmentsVector = 
00059     new G4PreCompoundFragmentVector(theFragmentsFactory->GetFragmentVector());
00060   g4pow = G4Pow::GetInstance();
00061   theParameters = G4PreCompoundParameters::GetAddress();
00062 }
00063 
00064 G4PreCompoundEmission::~G4PreCompoundEmission()
00065 {
00066   if (theFragmentsFactory) { delete theFragmentsFactory; }
00067   if (theFragmentsVector)  { delete theFragmentsVector; }
00068 }
00069 
00070 void G4PreCompoundEmission::SetDefaultModel()
00071 {
00072   if (theFragmentsFactory) { delete theFragmentsFactory; }
00073   theFragmentsFactory = new G4PreCompoundEmissionFactory();
00074   if (theFragmentsVector) 
00075     {
00076       theFragmentsVector->SetVector(theFragmentsFactory->GetFragmentVector());
00077     }
00078   else 
00079     {
00080       theFragmentsVector = 
00081         new G4PreCompoundFragmentVector(theFragmentsFactory->GetFragmentVector());
00082     }
00083   return;
00084 }
00085 
00086 void G4PreCompoundEmission::SetHETCModel()
00087 {
00088   if (theFragmentsFactory) delete theFragmentsFactory;
00089   theFragmentsFactory = new G4HETCEmissionFactory();
00090   if (theFragmentsVector) 
00091     {
00092       theFragmentsVector->SetVector(theFragmentsFactory->GetFragmentVector());
00093     }
00094   else 
00095     {
00096       theFragmentsVector = 
00097         new G4PreCompoundFragmentVector(theFragmentsFactory->GetFragmentVector());
00098     }
00099   return;
00100 }
00101 
00102 G4ReactionProduct * G4PreCompoundEmission::PerformEmission(G4Fragment & aFragment)
00103 {
00104   // Choose a Fragment for emission
00105   G4VPreCompoundFragment * thePreFragment = theFragmentsVector->ChooseFragment();
00106   if (thePreFragment == 0)
00107     {
00108       G4cout <<  "G4PreCompoundEmission::PerformEmission : I couldn't choose a fragment\n"
00109              << "while trying to de-excite\n" 
00110              << aFragment << G4endl;
00111       throw G4HadronicException(__FILE__, __LINE__, "");
00112     }
00113 
00114   //G4cout << "Chosen fragment: " << G4endl;
00115   //G4cout << *thePreFragment << G4endl;
00116 
00117   // Kinetic Energy of emitted fragment
00118   G4double kinEnergyOfEmittedFragment = thePreFragment->GetKineticEnergy(aFragment);
00119   //  if(kinEnergyOfEmittedFragment < MeV) {
00120   //  G4cout << "Chosen fragment: " << G4endl;
00121   //  G4cout << *thePreFragment << G4endl;
00122   //  G4cout << "Ekin= " << kinEnergyOfEmittedFragment << G4endl;
00123     // }
00124   if(kinEnergyOfEmittedFragment < 0.0) { kinEnergyOfEmittedFragment = 0.0; }
00125   
00126   // Calculate the fragment momentum (three vector)
00127   AngularDistribution(thePreFragment,aFragment,kinEnergyOfEmittedFragment);
00128   
00129   // Mass of emittef fragment
00130   G4double EmittedMass = thePreFragment->GetNuclearMass();
00131   // Now we can calculate the four momentum 
00132   // both options are valid and give the same result but 2nd one is faster
00133   G4LorentzVector Emitted4Momentum(theFinalMomentum,
00134                                    EmittedMass + kinEnergyOfEmittedFragment);
00135     
00136   // Perform Lorentz boost
00137   G4LorentzVector Rest4Momentum = aFragment.GetMomentum();
00138   Emitted4Momentum.boost(Rest4Momentum.boostVector());  
00139 
00140   // Set emitted fragment momentum
00141   thePreFragment->SetMomentum(Emitted4Momentum);        
00142 
00143   // NOW THE RESIDUAL NUCLEUS
00144   // ------------------------
00145 
00146   Rest4Momentum -= Emitted4Momentum;
00147     
00148   // Update nucleus parameters:
00149   // --------------------------
00150 
00151   // Z and A
00152   aFragment.SetZandA_asInt(thePreFragment->GetRestZ(),
00153                            thePreFragment->GetRestA());
00154     
00155   // Number of excitons
00156   aFragment.SetNumberOfParticles(aFragment.GetNumberOfParticles()-
00157                                  thePreFragment->GetA());
00158   // Number of charges
00159   aFragment.SetNumberOfCharged(aFragment.GetNumberOfCharged()-
00160                                thePreFragment->GetZ());
00161     
00162   // Update nucleus momentum 
00163   // A check on consistence of Z, A, and mass will be performed
00164   aFragment.SetMomentum(Rest4Momentum);
00165         
00166   // Create a G4ReactionProduct 
00167   G4ReactionProduct * MyRP = thePreFragment->GetReactionProduct();
00168 
00169   //  if(kinEnergyOfEmittedFragment < MeV) {
00170   //  G4cout << "G4PreCompoundEmission::Fragment emitted" << G4endl;
00171   //  G4cout << thePreFragment << G4endl;
00172     // }
00173   return MyRP;
00174 }
00175 
00176 void
00177 G4PreCompoundEmission::AngularDistribution(G4VPreCompoundFragment* thePreFragment,
00178                                            const G4Fragment& aFragment,
00179                                            G4double ekin) 
00180 {
00181   G4int p = aFragment.GetNumberOfParticles();
00182   G4int h = aFragment.GetNumberOfHoles();
00183   G4double U = aFragment.GetExcitationEnergy();
00184 
00185   // Emission particle separation energy
00186   G4double Bemission = thePreFragment->GetBindingEnergy();
00187 
00188   // Fermi energy
00189   G4double Ef = theParameters->GetFermiEnergy();
00190         
00191   //
00192   //  G4EvaporationLevelDensityParameter theLDP;
00193   //  G4double gg = (6.0/pi2)*aFragment.GetA()*
00194 
00195   G4double gg = (6.0/pi2)*aFragment.GetA_asInt()*theParameters->GetLevelDensity();
00196         
00197   // Average exciton energy relative to bottom of nuclear well
00198   G4double Eav = 2*p*(p+1)/((p+h)*gg);
00199         
00200   // Excitation energy relative to the Fermi Level
00201   G4double Uf = std::max(U - (p - h)*Ef , 0.0);
00202   //  G4double Uf = U - KineticEnergyOfEmittedFragment - Bemission;
00203 
00204   G4double w_num = rho(p+1, h, gg, Uf, Ef);
00205   G4double w_den = rho(p,   h, gg, Uf, Ef);
00206   if (w_num > 0.0 && w_den > 0.0)
00207     {
00208       Eav *= (w_num/w_den);
00209       Eav += - Uf/(p+h) + Ef;
00210     }
00211   else 
00212     {
00213       Eav = Ef;
00214     }
00215   
00216   // VI + JMQ 19/01/2010 update computation of the parameter an
00217   //
00218   G4double an = 0.0;
00219   G4double Eeff = ekin + Bemission + Ef;
00220   if(ekin > DBL_MIN && Eeff > DBL_MIN) {
00221 
00222     G4double zeta = std::max(1.0,9.3/std::sqrt(ekin/MeV));
00223   
00224     // This should be the projectile energy. If I would know which is 
00225     // the projectile (proton, neutron) I could remove the binding energy. 
00226     // But, what happens if INC precedes precompound? This approximation
00227     // seems to work well enough
00228     G4double ProjEnergy = aFragment.GetExcitationEnergy();
00229 
00230     an = 3*std::sqrt((ProjEnergy+Ef)*Eeff)/(zeta*Eav);
00231 
00232     G4int ne = aFragment.GetNumberOfExcitons() - 1;
00233     if ( ne > 1 ) { an /= (G4double)ne; }
00234                         
00235     // protection of exponent
00236     if ( an > 10. ) { an = 10.; }
00237   }
00238 
00239   // sample cosine of theta and not theta as in old versions  
00240   G4double random = G4UniformRand();
00241   G4double cost;
00242  
00243   if(an < 0.1) { cost = 1. - 2*random; }
00244   else {
00245     G4double exp2an = std::exp(-2*an);
00246     cost = 1. + std::log(1-random*(1-exp2an))/an;
00247     if(cost > 1.) { cost = 1.; }
00248     else if(cost < -1.) {cost = -1.; }
00249   }  
00250 
00251   G4double phi = CLHEP::twopi*G4UniformRand();
00252   
00253   // Calculate the momentum magnitude of emitted fragment       
00254   G4double pmag = std::sqrt(ekin*(ekin + 2.0*thePreFragment->GetNuclearMass()));
00255   
00256   G4double sint = std::sqrt((1.0-cost)*(1.0+cost));
00257 
00258   theFinalMomentum.set(pmag*std::cos(phi)*sint,pmag*std::sin(phi)*sint,pmag*cost);
00259 
00260   // theta is the angle wrt the incident direction
00261   G4ThreeVector theIncidentDirection = aFragment.GetMomentum().vect().unit();
00262   theFinalMomentum.rotateUz(theIncidentDirection);
00263 }
00264 
00265 G4double G4PreCompoundEmission::rho(G4int p, G4int h, G4double gg, 
00266                                     G4double E, G4double Ef) const
00267 {       
00268   // 25.02.2010 V.Ivanchenko added more protections
00269   G4double Aph   = (p*p + h*h + p - 3.0*h)/(4.0*gg);
00270   //  G4double alpha = (p*p + h*h)/(2.0*gg);
00271   
00272   if ( E - Aph < 0.0) { return 0.0; }
00273   
00274   G4double logConst =  (p+h)*std::log(gg) 
00275     - g4pow->logfactorial(p+h-1) - g4pow->logfactorial(p) - g4pow->logfactorial(h);
00276 
00277   // initialise values using j=0
00278 
00279   G4double t1=1;
00280   G4double t2=1;
00281   G4double logt3 = (p+h-1) * std::log(E-Aph) + logConst;
00282   const G4double logmax = 200.;
00283   if(logt3 > logmax) { logt3 = logmax; }
00284   G4double tot = std::exp( logt3 );
00285 
00286   // and now sum rest of terms
00287   // 25.02.2010 V.Ivanchenko change while to for loop and cleanup 
00288   G4double Eeff = E - Aph; 
00289   for(G4int j=1; j<=h; ++j) 
00290     {
00291       Eeff -= Ef;
00292       if(Eeff < 0.0) { break; }
00293       t1 *= -1.;
00294       t2 *= (G4double)(h+1-j)/(G4double)j;
00295       logt3 = (p+h-1) * std::log( Eeff) + logConst;
00296       if(logt3 > logmax) { logt3 = logmax; }
00297       tot += t1*t2*std::exp(logt3);
00298     }
00299         
00300   return tot;
00301 }

Generated on Mon May 27 17:49:24 2013 for Geant4 by  doxygen 1.4.7