source: trunk/source/processes/hadronic/models/low_energy/src/G4LEDeuteronInelastic.cc @ 1358

Last change on this file since 1358 was 819, checked in by garnier, 16 years ago

import all except CVS

File size: 4.6 KB
Line 
1//
2// ********************************************************************
3// * License and Disclaimer                                           *
4// *                                                                  *
5// * The  Geant4 software  is  copyright of the Copyright Holders  of *
6// * the Geant4 Collaboration.  It is provided  under  the terms  and *
7// * conditions of the Geant4 Software License,  included in the file *
8// * LICENSE and available at  http://cern.ch/geant4/license .  These *
9// * include a list of copyright holders.                             *
10// *                                                                  *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work  make  any representation or  warranty, express or implied, *
14// * regarding  this  software system or assume any liability for its *
15// * use.  Please see the license in the file  LICENSE  and URL above *
16// * for the full disclaimer and the limitation of liability.         *
17// *                                                                  *
18// * This  code  implementation is the result of  the  scientific and *
19// * technical work of the GEANT4 collaboration.                      *
20// * By using,  copying,  modifying or  distributing the software (or *
21// * any work based  on the software)  you  agree  to acknowledge its *
22// * use  in  resulting  scientific  publications,  and indicate your *
23// * acceptance of all terms of the Geant4 Software license.          *
24// ********************************************************************
25//
26//
27//
28 // Hadronic Process: Deuteron Inelastic Process
29 // J.L. Chuma, TRIUMF, 25-Feb-1997
30 // Last modified: 27-Mar-1997
31 // J.L. Chuma, 08-May-2001: Update original incident passed back in vec[0]
32 //                          from NuclearReaction
33 //
34#include "G4LEDeuteronInelastic.hh"
35#include "Randomize.hh"
36#include "G4Electron.hh"
37
38 G4HadFinalState *
39  G4LEDeuteronInelastic::ApplyYourself( const G4HadProjectile &aTrack,
40                                        G4Nucleus &targetNucleus )
41  { 
42    theParticleChange.Clear();
43    const G4HadProjectile *originalIncident = &aTrack;
44   
45    if( verboseLevel > 1 )
46    {
47      const G4Material *targetMaterial = aTrack.GetMaterial();
48      G4cout << "G4LEDeuteronInelastic::ApplyYourself called" << G4endl;
49      G4cout << "kinetic energy = " << originalIncident->GetKineticEnergy()/MeV << "MeV, ";
50      G4cout << "target material = " << targetMaterial->GetName() << ", ";
51    }
52   
53    // Work-around for lack of model above 100 MeV
54    if( originalIncident->GetKineticEnergy()/MeV > 100. ||
55        originalIncident->GetKineticEnergy() <= 0.1*MeV )
56    {
57      theParticleChange.SetStatusChange(isAlive);
58      theParticleChange.SetEnergyChange(aTrack.GetKineticEnergy());
59      theParticleChange.SetMomentumChange(aTrack.Get4Momentum().vect().unit()); 
60      return &theParticleChange;     
61    }
62
63    G4double A = targetNucleus.GetN();
64    G4double Z = targetNucleus.GetZ();
65    G4double theAtomicMass = targetNucleus.AtomicMass( A, Z );
66    G4double massVec[9];
67    massVec[0] = targetNucleus.AtomicMass( A+2.0, Z+1.0 );
68    massVec[1] = targetNucleus.AtomicMass( A+1.0, Z+1.0 );
69    massVec[2] = targetNucleus.AtomicMass( A+1.0, Z     );
70    massVec[3] = theAtomicMass;
71    massVec[4] = 0.;
72    if (A > 1.0 && A-1.0 > Z) 
73        massVec[4] = targetNucleus.AtomicMass( A-1.0, Z     );
74    massVec[5] = 0.;
75    if (A > 2.0 && Z > 1.0 && A-2.0 > Z-1.0) 
76        massVec[5] = targetNucleus.AtomicMass( A-2.0, Z-1.0 );
77    massVec[6] = 0.;
78    if (A > Z+1.0) 
79        massVec[6] = targetNucleus.AtomicMass( A,     Z+1.0 );
80    massVec[7] = massVec[3];
81    massVec[8] = 0.;
82    if (Z > 1.0) massVec[8] = targetNucleus.AtomicMass( A,Z-1.0 );
83   
84    G4FastVector<G4ReactionProduct,4> vec;  // vec will contain the secondary particles
85    G4int vecLen = 0;
86    vec.Initialize( 0 );
87   
88    theReactionDynamics.NuclearReaction( vec, vecLen, originalIncident,
89                                         targetNucleus, theAtomicMass, massVec );
90    //
91    G4double p = vec[0]->GetMomentum().mag();
92    theParticleChange.SetMomentumChange( vec[0]->GetMomentum() * (1.0/p)  );
93    theParticleChange.SetEnergyChange( vec[0]->GetKineticEnergy() );
94    delete vec[0];
95    //
96    G4DynamicParticle *pd;
97    for( G4int i=1; i<vecLen; ++i )
98    {
99      pd = new G4DynamicParticle();
100      pd->SetDefinition( vec[i]->GetDefinition() );
101      pd->SetMomentum( vec[i]->GetMomentum() );
102      theParticleChange.AddSecondary( pd );
103      delete vec[i];
104    }
105    return &theParticleChange;
106  }
107 
108 /* end of file */
109 
Note: See TracBrowser for help on using the repository browser.