G4LENDInelastic Class Reference

#include <G4LENDInelastic.hh>

Inheritance diagram for G4LENDInelastic:

G4LENDModel G4HadronicInteraction

Public Member Functions

 G4LENDInelastic (G4ParticleDefinition *pd)
 ~G4LENDInelastic ()
G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)

Detailed Description

Definition at line 44 of file G4LENDInelastic.hh.


Constructor & Destructor Documentation

G4LENDInelastic::G4LENDInelastic ( G4ParticleDefinition pd  )  [inline]

Definition at line 49 of file G4LENDInelastic.hh.

References G4LENDModel::create_used_target_map(), and G4LENDModel::proj.

00050      :G4LENDModel( "LENDInelastic" ) 
00051      { 
00052         proj = pd; 
00053        
00054 //        theModelName = "LENDInelastic for "; 
00055 //        theModelName += proj->GetParticleName(); 
00056         create_used_target_map();
00057      };

G4LENDInelastic::~G4LENDInelastic (  )  [inline]

Definition at line 59 of file G4LENDInelastic.hh.

00059 {;};


Member Function Documentation

G4HadFinalState * G4LENDInelastic::ApplyYourself ( const G4HadProjectile aTrack,
G4Nucleus aTargetNucleus 
) [virtual]

Reimplemented from G4LENDModel.

Definition at line 31 of file G4LENDInelastic.cc.

References G4HadFinalState::AddSecondary(), G4HadFinalState::Clear(), G4Gamma::Gamma(), G4HadProjectile::Get4Momentum(), G4Nucleus::GetA_asInt(), G4HadProjectile::GetKineticEnergy(), G4HadProjectile::GetMaterial(), G4LENDManager::GetNucleusEncoding(), G4GIDI_target::getOthersFinalState(), G4ParticleTable::GetParticleTable(), G4Material::GetTemperature(), G4Nucleus::GetZ_asInt(), G4LENDModel::lend_manager, G4Neutron::Neutron(), G4Proton::Proton(), G4DynamicParticle::SetDefinition(), G4DynamicParticle::SetMomentum(), G4HadFinalState::SetStatusChange(), stopAndKill, G4HadronicInteraction::theParticleChange, and G4LENDModel::usedTarget_map.

00032 {
00033 
00034    G4ThreeVector proj_p = aTrack.Get4Momentum().vect();
00035 
00036    G4double temp = aTrack.GetMaterial()->GetTemperature();
00037 
00038    //G4int iZ = int ( aTarg.GetZ() );
00039    //G4int iA = int ( aTarg.GetN() );
00040    //migrate to integer A and Z (GetN_asInt returns number of neutrons in the nucleus since this) 
00041    G4int iZ = aTarg.GetZ_asInt();
00042    G4int iA = aTarg.GetA_asInt();
00043    //G4cout << "target: Z = " << iZ << " N = " << iA << G4endl;
00044 
00045    G4double ke = aTrack.GetKineticEnergy();
00046    //G4cout << "projectile: KE = " << ke/MeV << " [MeV]" << G4endl;
00047 
00048    G4HadFinalState* theResult = &theParticleChange;
00049    theResult->Clear();
00050 
00051    G4GIDI_target* aTarget = usedTarget_map.find( lend_manager->GetNucleusEncoding( iZ , iA ) )->second->GetTarget();
00052    std::vector<G4GIDI_Product>* products = aTarget->getOthersFinalState( ke*MeV, temp, NULL, NULL );
00053    if ( products != NULL ) 
00054    {
00055 
00056       G4ThreeVector psum(0);
00057 
00058       int totN = 0;
00059       for ( G4int j = 0; j < int( products->size() ); j++ ) 
00060       {
00061 
00062          G4int jZ = (*products)[j].Z; 
00063          G4int jA = (*products)[j].A; 
00064 
00065          //G4cout << "ZA = " << 1000 * (*products)[j].Z + (*products)[j].A << "  EK = "
00066          //     << (*products)[j].kineticEnergy
00067          //     << " px  " <<  (*products)[j].px
00068          //     << " py  " <<  (*products)[j].py
00069          //     << " pz  " <<  (*products)[j].pz
00070          //     << G4endl;
00071 
00072          G4DynamicParticle* theSec = new G4DynamicParticle;
00073 
00074          if ( jA == 1 && jZ == 1 )
00075          {
00076             theSec->SetDefinition( G4Proton::Proton() );
00077             totN += 1;
00078          }
00079          else if ( jA == 1 && jZ == 0 )
00080          {
00081             theSec->SetDefinition( G4Neutron::Neutron() );
00082             totN += 1;
00083          } 
00084          else if ( jZ > 0 )
00085          {
00086             if ( jA != 0 )
00087             {
00088                theSec->SetDefinition( G4ParticleTable::GetParticleTable()->FindIon( jZ , jA , 0 , 0 ) );
00089                totN += jA;
00090             }
00091             else 
00092             {
00093                theSec->SetDefinition( G4ParticleTable::GetParticleTable()->FindIon( jZ , iA+1-totN , 0 , 0 ) );
00094             }
00095          } 
00096          else
00097          {
00098             theSec->SetDefinition( G4Gamma::Gamma() );
00099          } 
00100 
00101          G4ThreeVector p( (*products)[j].px*MeV , (*products)[j].py*MeV , (*products)[j].pz*MeV ); 
00102          psum += p; 
00103          if ( p.mag() == 0 ) p = proj_p - psum;
00104 
00105          theSec->SetMomentum( p );
00106 
00107          theResult->AddSecondary( theSec );
00108       } 
00109    }
00110    delete products;
00111 
00112    theResult->SetStatusChange( stopAndKill );
00113 
00114    return theResult; 
00115 
00116 }


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