G4PreCompoundModel.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 // by V. Lara
00029 //
00030 // Modified:
00031 // 01.04.2008 J.M.Quesada Several changes. Soft cut-off switched off. 
00032 // 01.05.2008 J.M.Quesada Protection against non-physical preeq. 
00033 //                        transitional regime has been set
00034 // 03.09.2008 J.M.Quesada for external choice of inverse cross section option
00035 // 06.09.2008 J.M.Quesada Also external choices have been added for:
00036 //                      - superimposed Coulomb barrier (useSICB=true) 
00037 //                      - "never go back"  hipothesis (useNGB=true) 
00038 //                      - soft cutoff from preeq. to equlibrium (useSCO=true)
00039 //                      - CEM transition probabilities (useCEMtr=true)  
00040 // 20.08.2010 V.Ivanchenko Cleanup of the code: 
00041 //                      - integer Z and A;
00042 //                      - emission and transition classes created at initialisation
00043 //                      - options are set at initialisation
00044 //                      - do not use copy-constructors for G4Fragment  
00045 // 03.01.2012 V.Ivanchenko Added pointer to G4ExcitationHandler to the 
00046 //                         constructor
00047 
00048 #include "G4PreCompoundModel.hh"
00049 #include "G4PhysicalConstants.hh"
00050 #include "G4SystemOfUnits.hh"
00051 #include "G4PreCompoundEmission.hh"
00052 #include "G4PreCompoundTransitions.hh"
00053 #include "G4GNASHTransitions.hh"
00054 #include "G4ParticleDefinition.hh"
00055 #include "G4Proton.hh"
00056 #include "G4Neutron.hh"
00057 
00058 #include "G4NucleiProperties.hh"
00059 #include "G4PreCompoundParameters.hh"
00060 #include "Randomize.hh"
00061 #include "G4DynamicParticle.hh"
00062 #include "G4ParticleTypes.hh"
00063 #include "G4ParticleTable.hh"
00064 #include "G4LorentzVector.hh"
00065 
00066 G4PreCompoundModel::G4PreCompoundModel(G4ExcitationHandler* ptr) 
00067   : G4VPreCompoundModel(ptr,"PRECO"), useHETCEmission(false), 
00068     useGNASHTransition(false), OPTxs(3), useSICB(false), 
00069     useNGB(false), useSCO(false), useCEMtr(true), maxZ(3), maxA(5) 
00070                                       //maxZ(9), maxA(17)
00071 {
00072   if(!ptr) { SetExcitationHandler(new G4ExcitationHandler()); }
00073   theParameters = G4PreCompoundParameters::GetAddress();
00074 
00075   theEmission = new G4PreCompoundEmission();
00076   if(useHETCEmission) { theEmission->SetHETCModel(); }
00077   else { theEmission->SetDefaultModel(); }
00078   theEmission->SetOPTxs(OPTxs);
00079   theEmission->UseSICB(useSICB);
00080 
00081   if(useGNASHTransition) { theTransition = new G4GNASHTransitions; }
00082   else { theTransition = new G4PreCompoundTransitions(); }
00083   theTransition->UseNGB(useNGB);
00084   theTransition->UseCEMtr(useCEMtr);
00085 
00086   proton = G4Proton::Proton();
00087   neutron = G4Neutron::Neutron();
00088 }
00089 
00090 G4PreCompoundModel::~G4PreCompoundModel() 
00091 {
00092   delete theEmission;
00093   delete theTransition;
00094   delete GetExcitationHandler();
00095 }
00096 
00097 void G4PreCompoundModel::ModelDescription(std::ostream& outFile) const
00098 {
00099         outFile << "The GEANT4 precompound model is considered as an extension of the\n"
00100                 <<      "hadron kinetic model. It gives a possibility to extend the low energy range\n"
00101                 <<      "of the hadron kinetic model for nucleon-nucleus inelastic collision and it \n"
00102                 <<      "provides a ”smooth” transition from kinetic stage of reaction described by the\n"
00103                 <<      "hadron kinetic model to the equilibrium stage of reaction described by the\n"
00104                 <<      "equilibrium deexcitation models.\n"
00105                 <<      "The initial information for calculation of pre-compound nuclear stage\n"
00106                 <<      "consists of the atomic mass number A, charge Z of residual nucleus, its\n"
00107                 <<      "four momentum P0 , excitation energy U and number of excitons n, which equals\n"
00108                 <<      "the sum of the number of particles p (from them p_Z are charged) and the number of\n"
00109                 <<      "holes h.\n"
00110                 <<      "At the preequilibrium stage of reaction, we follow the exciton model approach in ref. [1],\n"
00111                 <<      "taking into account the competition among all possible nuclear transitions\n"
00112                 <<      "with ∆n = +2, −2, 0 (which are defined by their associated transition probabilities) and\n"
00113                 <<      "the emission of neutrons, protons, deutrons, thritium and helium nuclei (also defined by\n"
00114                 <<      "their associated emission  probabilities according to exciton model)\n"
00115                 <<      "\n"
00116                 <<      "[1] K.K. Gudima, S.G. Mashnik, V.D. Toneev, Nucl. Phys. A401 329 (1983)\n"
00117                 << std::endl;
00118 }
00120 
00121 G4HadFinalState* G4PreCompoundModel::ApplyYourself(const G4HadProjectile & thePrimary,
00122                                                    G4Nucleus & theNucleus)
00123 {  
00124   const G4ParticleDefinition* primary = thePrimary.GetDefinition();
00125   if(primary != neutron && primary != proton) {
00126     std::ostringstream errOs;
00127     errOs << "BAD primary type in G4PreCompoundModel: " 
00128           << primary->GetParticleName() <<G4endl;
00129     throw G4HadronicException(__FILE__, __LINE__, errOs.str());
00130   }
00131 
00132   G4int Zp = 0;
00133   G4int Ap = 1;
00134   if(primary == proton) { Zp = 1; }
00135 
00136   G4int A = theNucleus.GetA_asInt();
00137   G4int Z = theNucleus.GetZ_asInt();
00138    
00139   //G4cout << "### G4PreCompoundModel::ApplyYourself: A= " << A << " Z= " << Z
00140   //     << " Ap= " << Ap << " Zp= " << Zp << G4endl; 
00141   // 4-Momentum
00142   G4LorentzVector p = thePrimary.Get4Momentum();
00143   G4double mass = G4NucleiProperties::GetNuclearMass(A, Z);
00144   p += G4LorentzVector(0.0,0.0,0.0,mass);
00145   //G4cout << "Primary 4-mom " << p << "  mass= " << mass << G4endl;
00146 
00147   // prepare fragment
00148   G4Fragment anInitialState(A + Ap, Z + Zp, p);
00149 
00150   // projectile and target nucleons
00151   // Add nucleon on which interaction happens
00152   //++Ap;
00153   //if(A*G4UniformRand() <= G4double(Z)) { Zp += 1; }
00154   anInitialState.SetNumberOfExcitedParticle(2, 1);
00155   anInitialState.SetNumberOfHoles(1,0);
00156   //  anInitialState.SetNumberOfExcitedParticle(Ap, Zp);
00157   // anInitialState.SetNumberOfHoles(Ap,Zp);
00158 
00159   anInitialState.SetCreationTime(thePrimary.GetGlobalTime());
00160   
00161   // call excitation handler
00162   G4ReactionProductVector * result = DeExcite(anInitialState);
00163 
00164   // fill particle change
00165   theResult.Clear();
00166   theResult.SetStatusChange(stopAndKill);
00167   for(G4ReactionProductVector::iterator i= result->begin(); i != result->end(); ++i)
00168     {
00169       G4DynamicParticle * aNew = 
00170         new G4DynamicParticle((*i)->GetDefinition(),
00171                               (*i)->GetTotalEnergy(),
00172                               (*i)->GetMomentum());
00173       delete (*i);
00174       theResult.AddSecondary(aNew);
00175     }
00176   delete result;
00177   
00178   //return the filled particle change
00179   return &theResult;
00180 }
00181 
00183 
00184 G4ReactionProductVector* G4PreCompoundModel::DeExcite(G4Fragment& aFragment)
00185 {
00186   G4ReactionProductVector * Result = new G4ReactionProductVector;
00187   G4double Eex = aFragment.GetExcitationEnergy();
00188   G4int Z = aFragment.GetZ_asInt(); 
00189   G4int A = aFragment.GetA_asInt(); 
00190 
00191   //G4cout << "### G4PreCompoundModel::DeExcite" << G4endl;
00192   //G4cout << aFragment << G4endl;
00193  
00194   // Perform Equilibrium Emission 
00195   if ((Z < maxZ && A < maxA) || Eex < MeV /*|| Eex > 3.*MeV*A*/) {
00196     PerformEquilibriumEmission(aFragment, Result);
00197     return Result;
00198   }
00199   
00200   // main loop  
00201   for (;;) {
00202     
00203     //fragment++;
00204     //G4cout<<"-------------------------------------------------------------------"<<G4endl;
00205     //G4cout<<"Fragment number .. "<<fragment<<G4endl;
00206     
00207     // Initialize fragment according with the nucleus parameters
00208     //G4cout << "### Loop over fragment" << G4endl;
00209     //G4cout << aFragment << G4endl;
00210 
00211     theEmission->Initialize(aFragment);
00212     
00213     G4double gg = (6.0/pi2)*aFragment.GetA_asInt()*theParameters->GetLevelDensity();
00214     
00215     G4int EquilibriumExcitonNumber = 
00216       G4lrint(std::sqrt(2*gg*aFragment.GetExcitationEnergy()));
00217     //   
00218     //    G4cout<<"Neq="<<EquilibriumExcitonNumber<<G4endl;
00219     //
00220     // J. M. Quesada (Jan. 08)  equilibrium hole number could be used as preeq.
00221     // evap. delimiter (IAEA report)
00222     
00223     // Loop for transitions, it is performed while there are preequilibrium transitions.
00224     G4bool ThereIsTransition = false;
00225     
00226     //        G4cout<<"----------------------------------------"<<G4endl;
00227     //        G4double NP=aFragment.GetNumberOfParticles();
00228     //        G4double NH=aFragment.GetNumberOfHoles();
00229     //        G4double NE=aFragment.GetNumberOfExcitons();
00230     //        G4cout<<" Ex. Energy="<<aFragment.GetExcitationEnergy()<<G4endl;
00231     //        G4cout<<"N. excitons="<<NE<<"  N. Part="<<NP<<"N. Holes ="<<NH<<G4endl;
00232     //G4int transition=0;
00233     do {
00234       //transition++;
00235       //G4cout<<"transition number .."<<transition<<G4endl;
00236       //G4cout<<" n ="<<aFragment.GetNumberOfExcitons()<<G4endl;
00237       G4bool go_ahead = false;
00238       // soft cutoff criterium as an "ad-hoc" solution to force increase in  evaporation  
00239       G4int test = aFragment.GetNumberOfExcitons();
00240       if (test <= EquilibriumExcitonNumber) { go_ahead=true; }
00241 
00242       //J. M. Quesada (Apr. 08): soft-cutoff switched off by default
00243       if (useSCO && !go_ahead)
00244         {
00245           G4double x = G4double(test)/G4double(EquilibriumExcitonNumber) - 1;
00246           if( G4UniformRand() < 1.0 -  std::exp(-x*x/0.32) ) { go_ahead = true; }
00247         } 
00248         
00249       // JMQ: WARNING:  CalculateProbability MUST be called prior to Get methods !! 
00250       // (O values would be returned otherwise)
00251       G4double TotalTransitionProbability = 
00252         theTransition->CalculateProbability(aFragment);
00253       G4double P1 = theTransition->GetTransitionProb1();
00254       G4double P2 = theTransition->GetTransitionProb2();
00255       G4double P3 = theTransition->GetTransitionProb3();
00256       //G4cout<<"#0 P1="<<P1<<" P2="<<P2<<"  P3="<<P3<<G4endl;
00257       
00258       //J.M. Quesada (May 2008) Physical criterium (lamdas)  PREVAILS over 
00259       //                        approximation (critical exciton number)
00260       //V.Ivanchenko (May 2011) added check on number of nucleons
00261       //                        to send a fragment to FermiBreakUp
00262       if(!go_ahead || P1 <= P2+P3 || 
00263          (aFragment.GetZ_asInt() < maxZ && aFragment.GetA_asInt() < maxA) )        
00264         {
00265           //G4cout<<"#4 EquilibriumEmission"<<G4endl; 
00266           PerformEquilibriumEmission(aFragment,Result);
00267           return Result;
00268         }
00269       else 
00270         {
00271           G4double TotalEmissionProbability = 
00272             theEmission->GetTotalProbability(aFragment);
00273           //
00274           //G4cout<<"#1 TotalEmissionProbability="<<TotalEmissionProbability<<" Nex= " 
00275           //    <<aFragment.GetNumberOfExcitons()<<G4endl;
00276           //
00277           // Check if number of excitons is greater than 0
00278           // else perform equilibrium emission
00279           if (aFragment.GetNumberOfExcitons() <= 0) 
00280             {
00281               PerformEquilibriumEmission(aFragment,Result);
00282               return Result;
00283             }
00284             
00285           //J.M.Quesada (May 08) this has already been done in order to decide  
00286           //                     what to do (preeq-eq) 
00287           // Sum of all probabilities
00288           G4double TotalProbability = TotalEmissionProbability 
00289             + TotalTransitionProbability;
00290             
00291           // Select subprocess
00292           if (TotalProbability*G4UniformRand() > TotalEmissionProbability) 
00293             {
00294               //G4cout<<"#2 Transition"<<G4endl; 
00295               // It will be transition to state with a new number of excitons
00296               ThereIsTransition = true;         
00297               // Perform the transition
00298               theTransition->PerformTransition(aFragment);
00299             } 
00300           else 
00301             {
00302               //G4cout<<"#3 Emission"<<G4endl; 
00303               // It will be fragment emission
00304               ThereIsTransition = false;
00305               Result->push_back(theEmission->PerformEmission(aFragment));
00306             }
00307         }
00308     } while (ThereIsTransition);   // end of do loop
00309   } // end of for (;;) loop
00310   return Result;
00311 }
00312 
00314 //       Initialisation
00316 
00317 void G4PreCompoundModel::UseHETCEmission() 
00318 { 
00319   useHETCEmission = true; 
00320   theEmission->SetHETCModel();
00321 }
00322 
00323 void G4PreCompoundModel::UseDefaultEmission() 
00324 { 
00325   useHETCEmission = false; 
00326   theEmission->SetDefaultModel();
00327 }
00328 
00329 void G4PreCompoundModel::UseGNASHTransition() { 
00330   useGNASHTransition = true; 
00331   delete theTransition;
00332   theTransition = new G4GNASHTransitions;
00333   theTransition->UseNGB(useNGB);
00334   theTransition->UseCEMtr(useCEMtr);
00335 }
00336 
00337 void G4PreCompoundModel::UseDefaultTransition() { 
00338   useGNASHTransition = false; 
00339   delete theTransition;
00340   theTransition = new G4PreCompoundTransitions();
00341   theTransition->UseNGB(useNGB);
00342   theTransition->UseCEMtr(useCEMtr);
00343 }
00344 
00345 void G4PreCompoundModel::SetOPTxs(G4int opt) 
00346 { 
00347   OPTxs = opt; 
00348   theEmission->SetOPTxs(OPTxs);
00349 }
00350 
00351 void G4PreCompoundModel::UseSICB() 
00352 { 
00353   useSICB = true; 
00354   theEmission->UseSICB(useSICB);
00355 }
00356 
00357 void G4PreCompoundModel::UseNGB()  
00358 { 
00359   useNGB = true; 
00360 }
00361 
00362 void G4PreCompoundModel::UseSCO()  
00363 { 
00364   useSCO = true; 
00365 }
00366 
00367 void G4PreCompoundModel::UseCEMtr() 
00368 { 
00369   useCEMtr = true; 
00370 }
00371 

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