source: trunk/source/processes/hadronic/stopping/include/G4AntiNeutronAnnihilationAtRest.hh @ 819

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

import all except CVS

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