Geant4  10.00.p03
G4VPreCompoundFragment.hh
Go to the documentation of this file.
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 // $Id: G4VPreCompoundFragment.hh 68028 2013-03-13 13:48:15Z gcosmo $
27 //
28 // J. M. Quesada (August 2008).
29 // Based on previous work by V. Lara
30 //
31 // Modif (03 September 2008) by J. M. Quesada for external choice of inverse
32 // cross section option
33 // JMQ (06 September 2008) Also external choice has been added for:
34 // - superimposed Coulomb barrier (if useSICB=true)
35 // 20.08.2010 V.Ivanchenko added int Z and A and cleanup; added
36 // G4ParticleDefinition to constructor,
37 // inline method to build G4ReactionProduct;
38 // remove string name
39 //
40 
41 #ifndef G4VPreCompoundFragment_h
42 #define G4VPreCompoundFragment_h 1
43 
44 #include "G4ios.hh"
45 #include <iomanip>
46 #include "G4ParticleDefinition.hh"
47 #include "G4IonTable.hh"
48 #include "G4Fragment.hh"
49 #include "G4VCoulombBarrier.hh"
50 #include "G4ReactionProduct.hh"
52 #include "G4Pow.hh"
53 
55 {
56 public:
57 
58  // ============================
59  // Constructors and destructor
60  // ============================
61 
63  G4VCoulombBarrier * aCoulombBarrier);
64 
65  virtual ~G4VPreCompoundFragment();
66 
67  // ==========
68  // operators
69  // ==========
70 
71  friend std::ostream&
72  operator<<(std::ostream&, const G4VPreCompoundFragment*);
73  friend std::ostream&
74  operator<<(std::ostream&, const G4VPreCompoundFragment&);
75 
76  // =====================
77  // Pure Virtual methods
78  // =====================
79 
80  // Initialization method
81  void Initialize(const G4Fragment & aFragment);
82 
83  // Methods for calculating the emission probability
84  // ------------------------------------------------
85 
86  // Calculates the total (integrated over kinetic energy) emission
87  // probability of a fragment
88  virtual G4double CalcEmissionProbability(const G4Fragment & aFragment) = 0;
89 
90  virtual G4double GetKineticEnergy(const G4Fragment & aFragment) = 0;
91 
92  inline G4ReactionProduct * GetReactionProduct() const;
93 
94  inline G4int GetA() const;
95 
96  inline G4int GetZ() const;
97 
98  inline G4int GetRestA() const;
99 
100  inline G4int GetRestZ() const;
101 
102  inline G4double ResidualA13() const;
103 
104  inline G4double GetCoulombBarrier() const;
105 
106  inline G4double GetBindingEnergy() const;
107 
108  inline G4double GetMaximalKineticEnergy() const;
109 
110  inline G4double GetEnergyThreshold() const;
111 
112  inline G4double GetEmissionProbability() const;
113 
114  inline G4double GetNuclearMass() const;
115 
116  inline G4double GetRestNuclearMass() const;
117 
118  inline G4double GetReducedMass() const;
119 
120  inline const G4LorentzVector& GetMomentum() const;
121 
122  inline void SetMomentum(const G4LorentzVector & value);
123 
124  inline const G4String GetName() const;
125 
126  //for inverse cross section choice
127  inline void SetOPTxs(G4int);
128  //for superimposed Coulomb Barrier for inverse cross sections
129  inline void UseSICB(G4bool);
130 
131 protected:
132 
133  inline G4bool IsItPossible(const G4Fragment & aFragment) const;
134 
135 private:
136 
137  // default constructor
139  // copy constructor
141  const G4VPreCompoundFragment&
142  operator= (const G4VPreCompoundFragment &right);
143  G4int operator==(const G4VPreCompoundFragment &right) const;
144  G4int operator!=(const G4VPreCompoundFragment &right) const;
145 
146  // =============
147  // Data members
148  // =============
149 
152 
157 
164 
166 
167 protected:
168 
171 
174 
175  //for inverse cross section choice
177  //for superimposed Coulomb Barrier for inverse cross sections
179 };
180 
182 
183 #endif
G4ReactionProduct * GetReactionProduct() const
G4double GetEnergyThreshold() const
Definition: G4Pow.hh:56
const G4LorentzVector & GetMomentum() const
G4double GetMaximalKineticEnergy() const
G4double ResidualA13() const
G4double GetReducedMass() const
const G4ParticleDefinition * particle
G4int operator!=(const G4VPreCompoundFragment &right) const
G4int GetA() const
G4bool IsItPossible(const G4Fragment &aFragment) const
int G4int
Definition: G4Types.hh:78
void Initialize(const G4Fragment &aFragment)
G4double GetBindingEnergy() const
void SetMomentum(const G4LorentzVector &value)
G4int operator==(const G4VPreCompoundFragment &right) const
bool G4bool
Definition: G4Types.hh:79
const G4String GetName() const
G4int GetRestZ() const
G4double GetCoulombBarrier() const
G4VCoulombBarrier * theCoulombBarrierPtr
G4double GetNuclearMass() const
G4double GetEmissionProbability() const
G4int GetRestA() const
G4PreCompoundParameters * theParameters
G4int GetZ() const
G4double GetRestNuclearMass() const
double G4double
Definition: G4Types.hh:76
virtual G4double GetKineticEnergy(const G4Fragment &aFragment)=0
virtual G4double CalcEmissionProbability(const G4Fragment &aFragment)=0
friend std::ostream & operator<<(std::ostream &, const G4VPreCompoundFragment *)
const G4VPreCompoundFragment & operator=(const G4VPreCompoundFragment &right)
CLHEP::HepLorentzVector G4LorentzVector