source: trunk/source/processes/hadronic/util/include/G4Nucleus.hh @ 1347

Last change on this file since 1347 was 1347, checked in by garnier, 13 years ago

geant4 tag 9.4

File size: 7.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 // original by H.P. Wellisch
29 // modified by J.L. Chuma, TRIUMF, 19-Nov-1996
30 // last modified: 27-Mar-1997
31 // Chr. Volcker, 10-Nov-1997: new methods and class variables.
32 // M.G. Pia, 2 Oct 1998: modified GetFermiMomentum (original design was
33 //                       the source of memory leaks)
34 // G.Folger, spring 2010:  add integer A/Z interface
35 
36#ifndef G4Nucleus_h
37#define G4Nucleus_h 1
38// Class Description
39// This class knows how to describe a nucleus;
40// to be used in your physics implementation (not physics list) in case you need this physics.
41// Class Description - End
42
43 
44#include "globals.hh"
45#include "G4ThreeVector.hh"
46#include "G4ParticleTypes.hh"
47#include "G4ReactionProduct.hh"
48#include "G4DynamicParticle.hh"
49#include "G4ReactionProductVector.hh"
50#include "Randomize.hh"
51 
52 class G4Nucleus 
53 {
54 public:
55   
56    G4Nucleus();
57   
58    G4Nucleus( const G4double A, const G4double Z );
59
60    G4Nucleus( const G4int A, const G4int Z );
61
62    G4Nucleus( const G4Material *aMaterial );
63   
64    ~G4Nucleus();
65   
66    inline G4Nucleus( const G4Nucleus &right )
67    { *this = right; }
68   
69    inline G4Nucleus & operator=( const G4Nucleus &right )
70     {
71       if( this != &right )
72       {
73         theA=right.theA; 
74         theZ=right.theZ; 
75         aEff=right.aEff; 
76         zEff=right.zEff; 
77         pnBlackTrackEnergy=right.pnBlackTrackEnergy; 
78         dtaBlackTrackEnergy=right.dtaBlackTrackEnergy;
79         pnBlackTrackEnergyfromAnnihilation =
80                      right.pnBlackTrackEnergyfromAnnihilation; 
81         dtaBlackTrackEnergyfromAnnihilation =
82                      right.dtaBlackTrackEnergyfromAnnihilation; 
83         theTemp = right.theTemp;
84         excitationEnergy = right.excitationEnergy;
85         momentum = right.momentum;
86         fermiMomentum = right.fermiMomentum;
87       }
88       return *this;
89     }
90   
91    inline G4bool operator==( const G4Nucleus &right ) const
92    { return ( this == (G4Nucleus *) &right ); }
93   
94    inline G4bool operator!=( const G4Nucleus &right ) const
95    { return ( this != (G4Nucleus *) &right ); }
96   
97    void ChooseParameters( const G4Material *aMaterial );
98
99    void SetParameters( const G4double A, const G4double Z );
100    void SetParameters( const G4int A, const G4int Z );
101   
102//deprecated Jan 2010, GF
103    inline G4double GetN() const
104    { return aEff; }
105   
106    inline G4double GetZ() const
107    { return zEff; }
108
109//to be replaced by new
110    inline G4int GetA_asInt() const
111    { return theA; }   
112   
113    inline G4int GetN_asInt() const
114    { return theA-theZ; }   
115   
116    inline G4int GetZ_asInt() const
117    { return theZ; }   
118//... \GF
119
120    G4DynamicParticle *ReturnTargetParticle() const;
121   
122    G4double AtomicMass( const G4double A, const G4double Z ) const;
123    G4double AtomicMass( const G4int A, const G4int Z ) const;
124   
125    G4double GetThermalPz( const G4double mass, const G4double temp ) const;
126   
127    G4ReactionProduct GetThermalNucleus(G4double aMass, G4double temp=-1) const;
128   
129    G4ReactionProduct GetBiasedThermalNucleus(G4double aMass, G4ThreeVector aVelocity, G4double temp=-1) const;
130
131    G4double Cinema( G4double kineticEnergy );
132   
133    G4double EvaporationEffects( G4double kineticEnergy );
134
135    G4double AnnihilationEvaporationEffects(G4double kineticEnergy, G4double ekOrg);
136   
137    inline G4double GetPNBlackTrackEnergy() const
138    { return pnBlackTrackEnergy; }
139   
140    inline G4double GetDTABlackTrackEnergy() const
141    { return dtaBlackTrackEnergy; }
142   
143    inline G4double GetAnnihilationPNBlackTrackEnergy() const
144    { return pnBlackTrackEnergyfromAnnihilation; }
145   
146    inline G4double GetAnnihilationDTABlackTrackEnergy() const
147    { return dtaBlackTrackEnergyfromAnnihilation; }
148   
149// ******************  methods introduced by ChV ***********************   
150   // return fermi momentum
151     G4ThreeVector GetFermiMomentum();
152
153/*
154  // return particle to be absorbed.
155     G4DynamicParticle* ReturnAbsorbingParticle(G4double weight);
156*/
157
158  //  final nucleus fragmentation. Return List of particles
159  // which should be used for further tracking.
160     G4ReactionProductVector* Fragmentate();
161     
162
163  // excitation Energy...
164     void AddExcitationEnergy(G4double anEnergy);
165 
166 
167  // momentum of absorbed Particles ..
168     void AddMomentum(const G4ThreeVector aMomentum);
169     
170  // return excitation Energy
171     G4double GetEnergyDeposit() {return excitationEnergy; }
172     
173
174
175// ****************************** end ChV ******************************
176
177
178 private:
179   
180    G4int    theA;
181    G4int    theZ;
182    G4double aEff;  // effective atomic weight
183    G4double zEff;  // effective atomic number
184   
185    G4double pnBlackTrackEnergy;  // the kinetic energy available for
186                                  // proton/neutron black track particles
187    G4double dtaBlackTrackEnergy; // the kinetic energy available for
188                                  // deuteron/triton/alpha particles
189    G4double pnBlackTrackEnergyfromAnnihilation;
190                     // kinetic energy available for proton/neutron black
191                     // track particles based on baryon annihilation
192    G4double dtaBlackTrackEnergyfromAnnihilation;
193                     // kinetic energy available for deuteron/triton/alpha
194                     // black track particles based on baryon annihilation
195
196
197// ************************** member variables by ChV *******************
198  // Excitation Energy leading to evaporation or deexcitation.
199     G4double  excitationEnergy;
200     
201  // Momentum, accumulated by absorbing Particles
202     G4ThreeVector momentum;
203     
204  // Fermi Gas model: at present, we assume constant nucleon density for all
205  // nuclei. The radius of a nucleon is taken to be 1 fm.
206  // see for example S.Fl"ugge, Encyclopedia of Physics, Vol XXXIX,
207  // Structure of Atomic Nuclei (Berlin-Gottingen-Heidelberg, 1957) page 426.
208
209  // maximum momentum possible from fermi gas model:
210     G4double fermiMomentum; 
211     G4double theTemp; // temperature
212// ****************************** end ChV ******************************
213
214 };
215 
216#endif
217 
Note: See TracBrowser for help on using the repository browser.