source: trunk/source/processes/hadronic/stopping/include/G4NeutronCaptureAtRest.hh @ 1055

Last change on this file since 1055 was 1055, checked in by garnier, 15 years ago

maj sur la beta de geant 4.9.3

File size: 4.0 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//   G4NeutronCaptureAtRest physics process
27//   Larry Felawka (TRIUMF), April 1998
28//---------------------------------------------------------------------
29
30#ifndef G4NeutronCaptureAtRest_h
31#define G4NeutronCaptureAtRest_h 1
32
33// Class Description:
34//
35// Process for capture of neutrons at rest.
36// To be used in your physics list in case you need this physics.
37
38 
39#include "globals.hh"
40#include "Randomize.hh"
41#include "G4VRestProcess.hh"
42#include "G4VParticleChange.hh"
43#include "G4ParticleDefinition.hh"
44#include "G4GHEKinematicsVector.hh"
45#include "G4HadronicProcessType.hh"
46
47class G4NeutronCaptureAtRest : public G4VRestProcess
48 
49{ 
50  private:
51  // hide assignment operator as private
52      G4NeutronCaptureAtRest& operator=(const G4NeutronCaptureAtRest &right);
53      G4NeutronCaptureAtRest(const G4NeutronCaptureAtRest& );
54   
55  public:
56 
57     G4NeutronCaptureAtRest(const G4String& processName ="NeutronCaptureAtRest", 
58                       G4ProcessType   aType = fHadronic );
59 
60    ~G4NeutronCaptureAtRest();
61
62     G4bool IsApplicable(const G4ParticleDefinition&);
63
64     void PreparePhysicsTable(const G4ParticleDefinition&);
65
66  // null physics table
67     void BuildPhysicsTable(const G4ParticleDefinition&);
68
69     G4double AtRestGetPhysicalInteractionLength(const G4Track&,
70                                                 G4ForceCondition*);
71
72  // zero mean lifetime
73     G4double GetMeanLifeTime(const G4Track& ,
74                              G4ForceCondition* ) {return 0.0;}
75
76     G4VParticleChange* AtRestDoIt(const G4Track&, const G4Step&); 
77
78  // return number of secondaries produced
79     G4int GetNumberOfSecondaries();
80
81  // pointer to array containg kinematics of secondaries
82     G4GHEKinematicsVector* GetSecondaryKinematics();
83
84  private:
85
86     void GenerateSecondaries();
87     void Normal( G4float* );
88     void NeutronCapture( G4int* );
89     G4double AtomAs( G4float, G4float );
90
91  private:
92
93// global time-of-flight of stopped hadron
94     G4float  globalTime;
95
96// atomic mass of target nucleus
97     G4float  targetAtomicMass;
98
99// charge of target nucleus
100     G4float  targetCharge;
101
102     G4GHEKinematicsVector* pv;
103     G4GHEKinematicsVector* eve;
104     G4GHEKinematicsVector* gkin;
105
106     G4int    ngkine;
107
108     G4int    ntot;
109     G4GHEKinematicsVector result;
110
111     G4float  massProton;
112     G4float  massNeutron;
113     G4float  massElectron;
114     G4float  massDeuteron;
115     G4float  massAlpha;
116
117     G4ParticleDefinition* pdefGamma;
118     G4ParticleDefinition* pdefNeutron;
119
120};
121
122#endif
123 
Note: See TracBrowser for help on using the repository browser.