G4LEAlphaInelastic.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 // Hadronic Process: Alpha Inelastic Process
00027 // J.L. Chuma, TRIUMF, 25-Feb-1997
00028 // J.L. Chuma, 08-May-2001: Update original incident passed back in vec[0]
00029 //                          from NuclearReaction
00030 
00031 #include <iostream>
00032 
00033 #include "G4LEAlphaInelastic.hh"
00034 #include "G4SystemOfUnits.hh"
00035 #include "Randomize.hh"
00036 #include "G4Electron.hh"
00037 
00038 G4LEAlphaInelastic::G4LEAlphaInelastic(const G4String& name)
00039  : G4InelasticInteraction(name)
00040 {
00041   SetMinEnergy(0.0*GeV);
00042   SetMaxEnergy(10.*TeV);
00043   G4cout << "WARNING: model G4LEAlphaInelastic is being deprecated and will\n"
00044          << "disappear in Geant4 version 10.0"  << G4endl;
00045 }
00046 
00047 
00048 void G4LEAlphaInelastic::ModelDescription(std::ostream& outFile) const
00049 {
00050   outFile << "G4LEAlphaInelastic is one of the Low Energy Parameterized\n"
00051           << "(LEP) models used to implement inelastic alpha scattering\n"
00052           << "from nuclei.  It is a re-engineered version of the GHEISHA\n"
00053           << "code of H. Fesefeldt. It divides the initial collision\n"
00054           << "products into backward- and forward-going clusters which are\n"
00055           << "then decayed into final state hadrons.  The model does not\n"
00056           << "conserve energy on an event-by-event basis.  It may be\n"
00057           << "applied to alphas with initial energies between 0 and 10\n"
00058           << "TeV.\n";
00059 }
00060 
00061 
00062 G4HadFinalState*
00063 G4LEAlphaInelastic::ApplyYourself(const G4HadProjectile& aTrack,
00064                                   G4Nucleus& targetNucleus)
00065 {
00066   theParticleChange.Clear();
00067   const G4HadProjectile* originalIncident = &aTrack;
00068 
00069   G4double A = targetNucleus.GetA_asInt();
00070   G4double Z = targetNucleus.GetZ_asInt();
00071         
00072   G4double kineticEnergy = aTrack.Get4Momentum().e()-aTrack.GetDefinition()->GetPDGMass();
00073   if (verboseLevel > 1) {
00074     const G4Material *targetMaterial = aTrack.GetMaterial();
00075     G4cout << "G4LEAlphaInelastic::ApplyYourself called" << G4endl;
00076     G4cout << "kinetc energy = " <<kineticEnergy/MeV << "MeV, ";
00077     G4cout << "target material = " << targetMaterial->GetName() << G4endl;
00078   }
00079     
00080   // Work-around for lack of model above 100 MeV
00081   if (kineticEnergy/MeV > 100. || kineticEnergy <= 0.1*MeV) {
00082     theParticleChange.SetStatusChange(isAlive);
00083     theParticleChange.SetEnergyChange(aTrack.GetKineticEnergy());
00084     theParticleChange.SetMomentumChange(aTrack.Get4Momentum().vect().unit()); 
00085     return &theParticleChange;      
00086   }
00087   G4double theAtomicMass = targetNucleus.AtomicMass( A, Z );
00088   G4double massVec[9];
00089   massVec[0] = targetNucleus.AtomicMass( A+4.0, Z+2.0 );
00090   massVec[1] = targetNucleus.AtomicMass( A+3.0, Z+2.0 );
00091   massVec[2] = targetNucleus.AtomicMass( A+3.0, Z+1.0 );
00092   massVec[3] = targetNucleus.AtomicMass( A+2.0, Z+1.0 );
00093   massVec[4] = targetNucleus.AtomicMass( A+1.0, Z+1.0 );
00094   massVec[5] = theAtomicMass;
00095   massVec[6] = targetNucleus.AtomicMass( A+2.0, Z+2.0 );
00096   massVec[7] = massVec[3];
00097   massVec[8] = targetNucleus.AtomicMass( A+2.0, Z     );
00098 
00099   G4FastVector<G4ReactionProduct,4> vec;  // vec will contain the secondary particles
00100   G4int vecLen = 0;
00101   vec.Initialize( 0 );
00102 
00103   theReactionDynamics.NuclearReaction(vec, vecLen, &aTrack,
00104                                       targetNucleus, theAtomicMass, massVec);
00105 
00106   G4double p = vec[0]->GetMomentum().mag();
00107   theParticleChange.SetMomentumChange( vec[0]->GetMomentum() *(1./p));
00108   theParticleChange.SetEnergyChange( vec[0]->GetKineticEnergy() );
00109   delete vec[0];
00110 
00111   if (vecLen <= 1) 
00112     {
00113       theParticleChange.SetStatusChange(isAlive);
00114       theParticleChange.SetEnergyChange(aTrack.GetKineticEnergy());
00115       theParticleChange.SetMomentumChange(aTrack.Get4Momentum().vect().unit()); 
00116       return &theParticleChange;      
00117     }
00118 
00119   G4DynamicParticle *pd;
00120   for (G4int i = 1; i < vecLen; ++i) {
00121     pd = new G4DynamicParticle();
00122     pd->SetDefinition( vec[i]->GetDefinition() );
00123     pd->SetMomentum( vec[i]->GetMomentum() );
00124     theParticleChange.AddSecondary( pd );
00125     delete vec[i];
00126   }
00127 
00128   if (isotopeProduction) DoIsotopeCounting(originalIncident, targetNucleus); 
00129   return &theParticleChange;
00130 }

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