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

Last change on this file since 1326 was 1315, checked in by garnier, 14 years ago

update geant4-09-04-beta-cand-01 interfaces-V09-03-09 vis-V09-03-08

File size: 7.4 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         theTemp = right.theTemp;
80         excitationEnergy = right.excitationEnergy;
81         momentum = right.momentum;
82         fermiMomentum = right.fermiMomentum;
83       }
84       return *this;
85     }
86   
87    inline G4bool operator==( const G4Nucleus &right ) const
88    { return ( this == (G4Nucleus *) &right ); }
89   
90    inline G4bool operator!=( const G4Nucleus &right ) const
91    { return ( this != (G4Nucleus *) &right ); }
92   
93    void ChooseParameters( const G4Material *aMaterial );
94
95    void SetParameters( const G4double A, const G4double Z );
96    void SetParameters( const G4int A, const G4int Z );
97   
98//deprecated Jan 2010, GF
99    inline G4double GetN() const
100    { return aEff; }
101   
102    inline G4double GetZ() const
103    { return zEff; }
104
105//to be replaced by new
106    inline G4int GetA_asInt() const
107    { return theA; }   
108   
109    inline G4int GetN_asInt() const
110    { return theA-theZ; }   
111   
112    inline G4int GetZ_asInt() const
113    { return theZ; }   
114//... \GF
115
116    G4DynamicParticle *ReturnTargetParticle() const;
117   
118    G4double AtomicMass( const G4double A, const G4double Z ) const;
119    G4double AtomicMass( const G4int A, const G4int Z ) const;
120   
121    G4double GetThermalPz( const G4double mass, const G4double temp ) const;
122   
123    G4ReactionProduct GetThermalNucleus(G4double aMass, G4double temp=-1) const;
124   
125    G4ReactionProduct GetBiasedThermalNucleus(G4double aMass, G4ThreeVector aVelocity, G4double temp=-1) const;
126
127    G4double Cinema( G4double kineticEnergy );
128   
129    G4double EvaporationEffects( G4double kineticEnergy );
130
131    G4double AnnihilationEvaporationEffects(G4double kineticEnergy, G4double ekOrg);
132   
133    inline G4double GetPNBlackTrackEnergy() const
134    { return pnBlackTrackEnergy; }
135   
136    inline G4double GetDTABlackTrackEnergy() const
137    { return dtaBlackTrackEnergy; }
138   
139    inline G4double GetAnnihilationPNBlackTrackEnergy() const
140    { return pnBlackTrackEnergyfromAnnihilation; }
141   
142    inline G4double GetAnnihilationDTABlackTrackEnergy() const
143    { return dtaBlackTrackEnergyfromAnnihilation; }
144   
145// ******************  methods introduced by ChV ***********************   
146   // return fermi momentum
147     G4ThreeVector GetFermiMomentum();
148
149/*
150  // return particle to be absorbed.
151     G4DynamicParticle* ReturnAbsorbingParticle(G4double weight);
152*/
153
154  //  final nucleus fragmentation. Return List of particles
155  // which should be used for further tracking.
156     G4ReactionProductVector* Fragmentate();
157     
158
159  // excitation Energy...
160     void AddExcitationEnergy(G4double anEnergy);
161 
162 
163  // momentum of absorbed Particles ..
164     void AddMomentum(const G4ThreeVector aMomentum);
165     
166  // return excitation Energy
167     G4double GetEnergyDeposit() {return excitationEnergy; }
168     
169
170
171// ****************************** end ChV ******************************
172
173
174 private:
175   
176    G4int    theA;
177    G4int    theZ;
178    G4double aEff;  // effective atomic weight
179    G4double zEff;  // effective atomic number
180   
181    G4double pnBlackTrackEnergy;  // the kinetic energy available for
182                                  // proton/neutron black track particles
183    G4double dtaBlackTrackEnergy; // the kinetic energy available for
184                                  // deuteron/triton/alpha particles
185    G4double pnBlackTrackEnergyfromAnnihilation;
186                     // kinetic energy available for proton/neutron black
187                     // track particles based on baryon annihilation
188    G4double dtaBlackTrackEnergyfromAnnihilation;
189                     // kinetic energy available for deuteron/triton/alpha
190                     // black track particles based on baryon annihilation
191
192
193// ************************** member variables by ChV *******************
194  // Excitation Energy leading to evaporation or deexcitation.
195     G4double  excitationEnergy;
196     
197  // Momentum, accumulated by absorbing Particles
198     G4ThreeVector momentum;
199     
200  // Fermi Gas model: at present, we assume constant nucleon density for all
201  // nuclei. The radius of a nucleon is taken to be 1 fm.
202  // see for example S.Fl"ugge, Encyclopedia of Physics, Vol XXXIX,
203  // Structure of Atomic Nuclei (Berlin-Gottingen-Heidelberg, 1957) page 426.
204
205  // maximum momentum possible from fermi gas model:
206     G4double fermiMomentum; 
207     G4double theTemp; // temperature
208// ****************************** end ChV ******************************
209
210 };
211 
212#endif
213 
Note: See TracBrowser for help on using the repository browser.