Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/hadronic/models/pre_equilibrium/exciton_model/include/G4PreCompoundEmission.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 // Hadronic Process: Nuclear Preequilibrium
 28 // 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 G4Pow and G4PreCompoundParameters pointers
 35 //                         use int Z and A and cleanup; 
 36 //                         inline methods moved from icc file to hh
 37 
 38 #ifndef G4PreCompoundEmission_h
 39 #define G4PreCompoundEmission_h 1
 40 
 41 #include "G4VPreCompoundFragment.hh"
 42 #include "G4ReactionProduct.hh"
 43 #include "G4Fragment.hh"
 44 #include "G4PreCompoundFragmentVector.hh"
 45 
 46 class G4VPreCompoundEmissionFactory;
 47 class G4Pow;
 48 class G4NuclearLevelData;
 49 
 50 class G4PreCompoundEmission
 51 {
 52 public:
 53 
 54   G4PreCompoundEmission();
 55 
 56   ~G4PreCompoundEmission();
 57 
 58   void SetDefaultModel();
 59 
 60   void SetHETCModel();
 61 
 62   G4ReactionProduct* PerformEmission(G4Fragment& aFragment);
 63 
 64   inline G4double GetTotalProbability(const G4Fragment& aFragment);
 65 
 66   inline void SetOPTxs(G4int);
 67 
 68   inline void UseSICB(G4bool);
 69 
 70   G4PreCompoundEmission(const G4PreCompoundEmission &right) = delete;
 71   const G4PreCompoundEmission& operator=
 72   (const G4PreCompoundEmission &right) = delete;
 73   G4bool operator==(const G4PreCompoundEmission &right) const = delete;
 74   G4bool operator!=(const G4PreCompoundEmission &right) const = delete;
 75 
 76 private:
 77 
 78   void AngularDistribution(G4VPreCompoundFragment * theFragment,
 79          const G4Fragment& aFragment,
 80          G4double kineticEnergy);
 81     
 82   G4double rho(G4int p, G4int h, G4double gg, 
 83                G4double E, G4double Ef) const;
 84 
 85   G4Pow* g4calc;
 86   G4NuclearLevelData* fNuclData;
 87 
 88   G4double fFermiEnergy;
 89 
 90   // A vector with the allowed emission fragments 
 91   G4PreCompoundFragmentVector* theFragmentsVector;
 92   G4VPreCompoundEmissionFactory* theFragmentsFactory;
 93 
 94   // Momentum of emitted fragment
 95   G4ThreeVector theFinalMomentum;
 96   G4bool fUseAngularGenerator;
 97 
 98   G4int fModelID;
 99 };
100 
101 inline G4double 
102 G4PreCompoundEmission::GetTotalProbability(const G4Fragment& aFragment) 
103 {
104   return theFragmentsVector->CalculateProbabilities(aFragment);
105 }
106 
107 inline void G4PreCompoundEmission::SetOPTxs(G4int opt)
108 {
109   theFragmentsVector->SetOPTxs(opt);
110 }
111 
112 inline void G4PreCompoundEmission::UseSICB(G4bool use)
113 {
114   theFragmentsVector->UseSICB(use);
115 }
116 
117 #endif
118