Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/hadronic/models/pre_equilibrium/exciton_model/include/G4VPreCompoundFragment.hh

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 // J. M. Quesada (August 2008).  
 28 // Based  on previous work by V. Lara
 29 //
 30 // Modif (03 September 2008) by J. M. Quesada for external choice of inverse 
 31 // cross section option 
 32 // JMQ (06 September 2008) Also external choice has been added for:
 33 //                      - superimposed Coulomb barrier (if useSICB=true) 
 34 // 20.08.2010 V.Ivanchenko added int Z and A and cleanup; added 
 35 //                        G4ParticleDefinition to constructor, 
 36 //                        inline method to build G4ReactionProduct; 
 37 //                        remove string name
 38 //                        
 39 
 40 #ifndef G4VPreCompoundFragment_h
 41 #define G4VPreCompoundFragment_h 1
 42 
 43 #include "G4ios.hh"
 44 #include <iomanip>
 45 #include "G4ParticleDefinition.hh"
 46 #include "G4IonTable.hh"
 47 #include "G4Fragment.hh"
 48 #include "G4ReactionProduct.hh"
 49 #include "G4Pow.hh"
 50 
 51 class G4NuclearLevelData;
 52 class G4DeexPrecoParameters;
 53 class G4VCoulombBarrier;
 54 class G4InterfaceToXS;
 55 
 56 class G4VPreCompoundFragment
 57 {
 58 public:  
 59 
 60   explicit G4VPreCompoundFragment(const G4ParticleDefinition*,
 61           G4VCoulombBarrier*);
 62   
 63   virtual ~G4VPreCompoundFragment();
 64   
 65   friend std::ostream& 
 66   operator<<(std::ostream&, const G4VPreCompoundFragment*);
 67   friend std::ostream& 
 68   operator<<(std::ostream&, const G4VPreCompoundFragment&);
 69   
 70   // =====================
 71   // Pure Virtual methods
 72   // =====================
 73   
 74   // Run time initialization method
 75   G4bool Initialize(const G4Fragment& aFragment);
 76     
 77   // Methods for calculating the emission probability
 78   // ------------------------------------------------
 79   
 80   // Calculates the total (integrated over kinetic energy) emission
 81   // probability of a fragment
 82   virtual G4double CalcEmissionProbability(const G4Fragment&) = 0;
 83   
 84   // sample kinetic energy of emitted fragment
 85   virtual G4double SampleKineticEnergy(const G4Fragment&) = 0;
 86 
 87   inline G4ReactionProduct* GetReactionProduct() const;   
 88   
 89   G4int GetA() const { return theA; }
 90   
 91   G4int GetZ() const { return theZ; }
 92   
 93   G4int GetRestA() const { return theResA; }
 94   
 95   G4int GetRestZ() const { return theResZ; }
 96 
 97   G4double GetBindingEnergy() const { return theBindingEnergy; }
 98   
 99   G4double GetEnergyThreshold() const
100   {
101     return theMaxKinEnergy - theCoulombBarrier;
102   }
103 
104   G4double GetEmissionProbability() const { return theEmissionProbability; }
105     
106   G4double GetNuclearMass() const { return theMass; }
107   
108   G4double GetRestNuclearMass() const { return theResMass; }
109 
110   const G4LorentzVector& GetMomentum() const { return theMomentum; }
111   
112   void SetMomentum(const G4LorentzVector& lv) { theMomentum = lv; }
113   
114   //for inverse cross section choice
115   void SetOPTxs(G4int) {}
116   //for superimposed Coulomb Barrier for inverse cross sections
117   void UseSICB(G4bool use) { useSICB = use; } 
118 
119   G4VPreCompoundFragment(const G4VPreCompoundFragment &right) = delete;
120   const G4VPreCompoundFragment& 
121   operator= (const G4VPreCompoundFragment &right) = delete;  
122   G4bool operator==(const G4VPreCompoundFragment &right) const = delete;
123   G4bool operator!=(const G4VPreCompoundFragment &right) const = delete;
124 
125 protected:
126 
127   virtual G4double GetAlpha() const = 0;
128 
129   virtual G4double GetBeta() const { return -theCoulombBarrier; }
130 
131   G4NuclearLevelData* fNucData;
132   G4DeexPrecoParameters* theParameters;
133   G4Pow* g4calc;
134   G4InterfaceToXS* fXSection{nullptr};
135 
136   G4int theA;
137   G4int theZ;
138   G4int theResA{0};
139   G4int theResZ{0};
140   G4int theFragA{0};
141   G4int theFragZ{0};
142   //for inverse cross section choice
143   G4int OPTxs;
144   G4int index{0};
145 
146   G4double theResA13{0.0};
147   G4double theBindingEnergy{0.0};
148   G4double theMinKinEnergy{0.0};
149   G4double theMaxKinEnergy{0.0};
150   G4double theResMass{0.0};
151   G4double theReducedMass{0.0};
152   G4double theMass;
153 
154   G4double theEmissionProbability{0.0};
155   G4double theCoulombBarrier{0.0};
156 
157   //for superimposed Coulomb Barrier for inverse cross sections
158   G4bool useSICB{true};
159 
160 private:
161 
162   const G4ParticleDefinition* particle;
163   G4VCoulombBarrier* theCoulombBarrierPtr;
164   G4LorentzVector theMomentum{0., 0., 0., 0.};
165 };
166 
167 inline G4ReactionProduct* G4VPreCompoundFragment::GetReactionProduct() const
168 {
169   G4ReactionProduct* theReactionProduct = new G4ReactionProduct(particle);
170   theReactionProduct->SetMomentum(GetMomentum().vect());
171   theReactionProduct->SetTotalEnergy(GetMomentum().e());
172   return theReactionProduct;
173 }
174 
175 #endif
176