G4INCLClusterDecay.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 // INCL++ intra-nuclear cascade model
00027 // Pekka Kaitaniemi, CEA and Helsinki Institute of Physics
00028 // Davide Mancusi, CEA
00029 // Alain Boudard, CEA
00030 // Sylvie Leray, CEA
00031 // Joseph Cugnon, University of Liege
00032 //
00033 #define INCLXX_IN_GEANT4_MODE 1
00034 
00035 #include "globals.hh"
00036 
00044 #include "G4INCLClusterDecay.hh"
00045 #include "G4INCLParticleTable.hh"
00046 #include "G4INCLKinematicsUtils.hh"
00047 #include "G4INCLRandom.hh"
00048 // #include <cassert>
00049 #include <algorithm>
00050 
00051 namespace G4INCL {
00052 
00053   ParticleList ClusterDecay::decay(Cluster * const c) {
00054     ParticleList decayProducts;
00055     recursiveDecay(c, &decayProducts);
00056 
00057     // Correctly update the particle type
00058     if(c->getA()==1) {
00059 // assert(c->getZ()==1 || c->getZ()==0);
00060       if(c->getZ()==1)
00061         c->setType(Proton);
00062       else
00063         c->setType(Neutron);
00064       c->setTableMass();
00065     }
00066 
00067     return decayProducts;
00068   }
00069 
00070   void ClusterDecay::recursiveDecay(Cluster * const c, ParticleList *decayProducts) {
00071     const G4int Z = c->getZ();
00072     const G4int A = c->getA();
00073 // assert(c->getExcitationEnergy()>-1.e-5);
00074     if(c->getExcitationEnergy()<0.)
00075       c->setExcitationEnergy(0.);
00076 
00077     if(Z<ParticleTable::clusterTableZSize && A<ParticleTable::clusterTableASize) {
00078       ParticleTable::ClusterDecayType theDecayMode = ParticleTable::clusterDecayMode[Z][A];
00079 
00080       switch(theDecayMode) {
00081         default:
00082           ERROR("Unrecognized cluster-decay mode: " << theDecayMode << std::endl
00083               << c->print());
00084         case ParticleTable::StableCluster:
00085           // For stable clusters, just return
00086           return;
00087           break;
00088         case ParticleTable::ProtonDecay:
00089         case ParticleTable::NeutronDecay:
00090         case ParticleTable::AlphaDecay:
00091           // Two-body decays
00092           twoBodyDecay(c, theDecayMode, decayProducts);
00093           break;
00094         case ParticleTable::TwoProtonDecay:
00095         case ParticleTable::TwoNeutronDecay:
00096           // Three-body decays
00097           threeBodyDecay(c, theDecayMode, decayProducts);
00098           break;
00099         case ParticleTable::ProtonUnbound:
00100         case ParticleTable::NeutronUnbound:
00101           // Phase-space decays
00102           phaseSpaceDecay(c, theDecayMode, decayProducts);
00103           break;
00104       }
00105 
00106       // Calls itself recursively in case the produced remnant is still unstable.
00107       // Sneaky, isn't it.
00108       recursiveDecay(c,decayProducts);
00109 
00110     } else {
00111       // The cluster is too large for our decay-mode table. Decompose it only
00112       // if Z==0 || Z==A.
00113       DEBUG("Cluster is outside the decay-mode table." << c->print() << std::endl);
00114       if(Z==A) {
00115         DEBUG("Z==A, will decompose it in free protons." << std::endl);
00116         phaseSpaceDecay(c, ParticleTable::ProtonUnbound, decayProducts);
00117       } else if(Z==0) {
00118         DEBUG("Z==0, will decompose it in free neutrons." << std::endl);
00119         phaseSpaceDecay(c, ParticleTable::NeutronUnbound, decayProducts);
00120       }
00121     }
00122   }
00123 
00124   void ClusterDecay::twoBodyDecay(Cluster * const c, ParticleTable::ClusterDecayType theDecayMode, ParticleList *decayProducts) {
00125     Particle *decayParticle = 0;
00126     const ThreeVector mom(0.0, 0.0, 0.0);
00127     const ThreeVector pos = c->getPosition();
00128 
00129     // Create the emitted particle
00130     switch(theDecayMode) {
00131       case ParticleTable::ProtonDecay:
00132         decayParticle = new Particle(Proton, mom, pos);
00133         break;
00134       case ParticleTable::NeutronDecay:
00135         decayParticle = new Particle(Neutron, mom, pos);
00136         break;
00137       case ParticleTable::AlphaDecay:
00138         decayParticle = new Cluster(2,4);
00139         break;
00140       default:
00141         ERROR("Unrecognized cluster-decay mode in two-body decay: " << theDecayMode << std::endl
00142             << c->print());
00143         return;
00144     }
00145     decayParticle->makeParticipant();
00146     decayParticle->setNumberOfDecays(1);
00147     decayParticle->setPosition(c->getPosition());
00148     decayParticle->setEmissionTime(c->getEmissionTime());
00149     decayParticle->setTableMass();
00150 
00151     // Save some variables of the mother cluster
00152     G4double motherMass = c->getMass();
00153     const ThreeVector velocity = -c->boostVector();
00154 
00155     // Characteristics of the daughter particle
00156     const G4int daughterZ = c->getZ() - decayParticle->getZ();
00157     const G4int daughterA = c->getA() - decayParticle->getA();
00158     const G4double daughterMass = ParticleTable::getTableMass(daughterA,daughterZ);
00159 
00160     // The mother cluster becomes the daughter
00161     c->setZ(daughterZ);
00162     c->setA(daughterA);
00163     c->setMass(daughterMass);
00164     c->setExcitationEnergy(0.);
00165 
00166     // Decay kinematics in the mother rest frame
00167     const G4double decayMass = decayParticle->getMass();
00168 // assert(motherMass-daughterMass-decayMass>-1.e-5); // Q-value should be >0
00169     G4double pCM = 0.;
00170     if(motherMass-daughterMass-decayMass>0.)
00171       pCM = KinematicsUtils::momentumInCM(motherMass, daughterMass, decayMass);
00172     const ThreeVector momentum = Random::normVector(pCM);
00173     c->setMomentum(momentum);
00174     c->adjustEnergyFromMomentum();
00175     decayParticle->setMomentum(-momentum);
00176     decayParticle->adjustEnergyFromMomentum();
00177 
00178     // Boost to the lab frame
00179     decayParticle->boost(velocity);
00180     c->boost(velocity);
00181 
00182     // Add the decay particle to the list of decay products
00183     decayProducts->push_back(decayParticle);
00184   }
00185 
00186   void ClusterDecay::threeBodyDecay(Cluster * const c, ParticleTable::ClusterDecayType theDecayMode, ParticleList *decayProducts) {
00187     Particle *decayParticle1 = 0;
00188     Particle *decayParticle2 = 0;
00189     const ThreeVector mom(0.0, 0.0, 0.0);
00190     const ThreeVector pos = c->getPosition();
00191 
00192     // Create the emitted particles
00193     switch(theDecayMode) {
00194       case ParticleTable::TwoProtonDecay:
00195         decayParticle1 = new Particle(Proton, mom, pos);
00196         decayParticle2 = new Particle(Proton, mom, pos);
00197         break;
00198       case ParticleTable::TwoNeutronDecay:
00199         decayParticle1 = new Particle(Neutron, mom, pos);
00200         decayParticle2 = new Particle(Neutron, mom, pos);
00201         break;
00202       default:
00203         ERROR("Unrecognized cluster-decay mode in three-body decay: " << theDecayMode << std::endl
00204             << c->print());
00205         return;
00206     }
00207     decayParticle1->makeParticipant();
00208     decayParticle2->makeParticipant();
00209     decayParticle1->setNumberOfDecays(1);
00210     decayParticle2->setNumberOfDecays(1);
00211     decayParticle1->setTableMass();
00212     decayParticle2->setTableMass();
00213 
00214     // Save some variables of the mother cluster
00215     const G4double motherMass = c->getMass();
00216     const ThreeVector velocity = -c->boostVector();
00217 
00218     // Masses and charges of the daughter particle and of the decay products
00219     const G4int decayZ1 = decayParticle1->getZ();
00220     const G4int decayA1 = decayParticle1->getA();
00221     const G4int decayZ2 = decayParticle2->getZ();
00222     const G4int decayA2 = decayParticle2->getA();
00223     const G4int decayZ = decayZ1 + decayZ2;
00224     const G4int decayA = decayA1 + decayA2;
00225     const G4int daughterZ = c->getZ() - decayZ;
00226     const G4int daughterA = c->getA() - decayA;
00227     const G4double decayMass1 = decayParticle1->getMass();
00228     const G4double decayMass2 = decayParticle2->getMass();
00229     const G4double daughterMass = ParticleTable::getTableMass(daughterA,daughterZ);
00230 
00231     // Q-values
00232     G4double qValue = motherMass - daughterMass - decayMass1 - decayMass2;
00233 // assert(qValue > -1e-5); // Q-value should be >0
00234     if(qValue<0.)
00235       qValue=0.;
00236     const G4double qValueB = qValue * Random::shoot();
00237 
00238     // The decay particles behave as if they had more mass until the second
00239     // decay
00240     const G4double decayMass = decayMass1 + decayMass2 + qValueB;
00241 
00242     /* Stage A: mother --> daughter + (decay1+decay2) */
00243 
00244     // The mother cluster becomes the daughter
00245     c->setZ(daughterZ);
00246     c->setA(daughterA);
00247     c->setMass(daughterMass);
00248     c->setExcitationEnergy(0.);
00249 
00250     // Decay kinematics in the mother rest frame
00251     const G4double pCMA = KinematicsUtils::momentumInCM(motherMass, daughterMass, decayMass);
00252     const ThreeVector momentumA = Random::normVector(pCMA);
00253     c->setMomentum(momentumA);
00254     c->adjustEnergyFromMomentum();
00255     const ThreeVector decayBoostVector = momentumA/std::sqrt(decayMass*decayMass + momentumA.mag2());
00256 
00257     /* Stage B: (decay1+decay2) --> decay1 + decay2 */
00258 
00259     // Decay kinematics in the (decay1+decay2) rest frame
00260     const G4double pCMB = KinematicsUtils::momentumInCM(decayMass, decayMass1, decayMass2);
00261     const ThreeVector momentumB = Random::normVector(pCMB);
00262     decayParticle1->setMomentum(momentumB);
00263     decayParticle2->setMomentum(-momentumB);
00264     decayParticle1->adjustEnergyFromMomentum();
00265     decayParticle2->adjustEnergyFromMomentum();
00266 
00267     // Boost decay1 and decay2 to the Stage-A decay frame
00268     decayParticle1->boost(decayBoostVector);
00269     decayParticle2->boost(decayBoostVector);
00270 
00271     // Boost all particles to the lab frame
00272     decayParticle1->boost(velocity);
00273     decayParticle2->boost(velocity);
00274     c->boost(velocity);
00275 
00276     // Add the decay particles to the list of decay products
00277     decayProducts->push_back(decayParticle1);
00278     decayProducts->push_back(decayParticle2);
00279   }
00280 
00281   void ClusterDecay::phaseSpaceDecay(Cluster * const c, ParticleTable::ClusterDecayType theDecayMode, ParticleList *decayProducts) {
00282     const G4int theA = c->getA();
00283     const G4int theZ = c->getZ();
00284     const ThreeVector mom(0.0, 0.0, 0.0);
00285     const ThreeVector pos = c->getPosition();
00286 
00287     G4int theZStep;
00288     ParticleType theEjectileType;
00289     switch(theDecayMode) {
00290       case ParticleTable::ProtonUnbound:
00291         theZStep = 1;
00292         theEjectileType = Proton;
00293         break;
00294       case ParticleTable::NeutronUnbound:
00295         theZStep = 0;
00296         theEjectileType = Neutron;
00297         break;
00298       default:
00299         ERROR("Unrecognized cluster-decay mode in phase-space decay: " << theDecayMode << std::endl
00300             << c->print());
00301         return;
00302     }
00303 
00304     // Find the daughter cluster (first cluster which is not
00305     // proton/neutron-unbound, in the sense of the table)
00306     G4int finalDaughterZ, finalDaughterA;
00307     if(theZ<ParticleTable::clusterTableZSize && theA<ParticleTable::clusterTableASize) {
00308       finalDaughterZ=theZ;
00309       finalDaughterA=theA;
00310       while(ParticleTable::clusterDecayMode[finalDaughterZ][finalDaughterA]==theDecayMode) {
00311         finalDaughterA--;
00312         finalDaughterZ -= theZStep;
00313       }
00314     } else {
00315       finalDaughterA = 1;
00316       if(theDecayMode==ParticleTable::ProtonUnbound)
00317         finalDaughterZ = 1;
00318       else
00319         finalDaughterZ = 0;
00320     }
00321 // assert(finalDaughterZ<=theZ && finalDaughterA<theA && finalDaughterA>0 && finalDaughterZ>=0);
00322     const G4double finalDaughterMass = ParticleTable::getTableMass(finalDaughterA, finalDaughterZ);
00323 
00324     // Compute the available decay energy
00325     const G4int nSplits = theA-finalDaughterA;
00326     const G4double ejectileMass = ParticleTable::getTableMass(1, theZStep);
00327     // c->getMass() can possibly contain some excitation energy, too
00328     G4double availableEnergy = c->getMass() - finalDaughterMass - nSplits*ejectileMass;
00329 // assert(availableEnergy>-1.e-5);
00330     if(availableEnergy<0.)
00331       availableEnergy=0.;
00332 
00333     // Compute an estimate of the maximum event weight
00334     G4double maximumWeight = 1.;
00335     G4double eMax = finalDaughterMass + availableEnergy;
00336     G4double eMin = finalDaughterMass - ejectileMass;
00337     for(G4int iSplit=0; iSplit<nSplits; ++iSplit) {
00338       eMax += ejectileMass;
00339       eMin += ejectileMass;
00340       maximumWeight *= KinematicsUtils::momentumInCM(eMax, eMin, ejectileMass);
00341     }
00342 
00343     // Sample decays until the weight cutoff is satisfied
00344     G4double weight;
00345     std::vector<G4double> theCMMomenta;
00346     std::vector<G4double> invariantMasses;
00347     G4int nTries=0;
00348     /* Maximum number of trials dependent on nSplits. 50 trials seems to be
00349      * sufficient for small nSplits. For nSplits>=5, maximumWeight is a gross
00350      * overestimate of the actual maximum weight, leading to unreasonably high
00351      * rejection rates. For these cases, we set nSplits=1000, although the sane
00352      * thing to do would be to improve the importance sampling (maybe by
00353      * parametrising maximumWeight?).
00354      */
00355     G4int maxTries;
00356     if(nSplits<5)
00357       maxTries=50;
00358     else
00359       maxTries=1000;
00360     do {
00361       if(nTries++>maxTries) {
00362         WARN("Phase-space decay exceeded the maximum number of rejections (" << maxTries
00363             << "). Z=" << theZ << ", A=" << theA << ", E*=" << c->getExcitationEnergy()
00364             << ", availableEnergy=" << availableEnergy
00365             << ", nSplits=" << nSplits
00366             << std::endl);
00367         break;
00368       }
00369 
00370       // Construct a sorted vector of random numbers
00371       std::vector<G4double> randomNumbers;
00372       for(G4int iSplit=0; iSplit<nSplits-1; ++iSplit)
00373         randomNumbers.push_back(Random::shoot0());
00374       std::sort(randomNumbers.begin(), randomNumbers.end());
00375 
00376       // Divide the available decay energy in the right number of steps
00377       invariantMasses.clear();
00378       invariantMasses.push_back(finalDaughterMass);
00379       for(G4int iSplit=0; iSplit<nSplits-1; ++iSplit)
00380         invariantMasses.push_back(finalDaughterMass + (iSplit+1)*ejectileMass + randomNumbers.at(iSplit)*availableEnergy);
00381       invariantMasses.push_back(c->getMass());
00382 
00383       weight = 1.;
00384       theCMMomenta.clear();
00385       for(G4int iSplit=0; iSplit<nSplits; ++iSplit) {
00386         G4double motherMass = invariantMasses.at(nSplits-iSplit);
00387         const G4double daughterMass = invariantMasses.at(nSplits-iSplit-1);
00388 // assert(motherMass-daughterMass-ejectileMass>-1.e-5);
00389         G4double pCM = 0.;
00390         if(motherMass-daughterMass-ejectileMass>0.)
00391           pCM = KinematicsUtils::momentumInCM(motherMass, daughterMass, ejectileMass);
00392         theCMMomenta.push_back(pCM);
00393         weight *= pCM;
00394       }
00395     } while(maximumWeight*Random::shoot()>weight);
00396 
00397     for(G4int iSplit=0; iSplit<nSplits; ++iSplit) {
00398       ThreeVector const velocity = -c->boostVector();
00399 
00400 #if !defined(NDEBUG) && !defined(INCLXX_IN_GEANT4_MODE)
00401       const G4double motherMass = c->getMass();
00402 #endif
00403       c->setA(c->getA() - 1);
00404       c->setZ(c->getZ() - theZStep);
00405       c->setMass(invariantMasses.at(nSplits-iSplit-1));
00406 
00407       Particle *ejectile = new Particle(theEjectileType, mom, pos);
00408       ejectile->setTableMass();
00409 
00410 // assert(motherMass-c->getMass()-ejectileMass>-1.e-5);
00411       ThreeVector momentum;
00412       momentum = Random::normVector(theCMMomenta.at(iSplit));
00413       c->setMomentum(momentum);
00414       c->adjustEnergyFromMomentum();
00415       ejectile->setMomentum(-momentum);
00416       ejectile->adjustEnergyFromMomentum();
00417 
00418       // Boost to the lab frame
00419       c->boost(velocity);
00420       ejectile->boost(velocity);
00421 
00422       // Add the decay particle to the list of decay products
00423       decayProducts->push_back(ejectile);
00424     }
00425 // assert(std::abs(c->getTableMass()-c->getMass())<1.e-3);
00426     c->setExcitationEnergy(0.);
00427   }
00428 }
00429 

Generated on Mon May 27 17:48:34 2013 for Geant4 by  doxygen 1.4.7