Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/hadronic/models/pre_equilibrium/exciton_model/src/G4VPreCompoundFragment.cc

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).  Based  on previous work by V. Lara
 28 //
 29 // Modified:
 30 // 20.08.2010 V.Ivanchenko added G4Pow and G4PreCompoundParameters pointers
 31 //                         use int Z and A and cleanup
 32 
 33 #include "G4VPreCompoundFragment.hh"
 34 #include "G4SystemOfUnits.hh"
 35 #include "G4NucleiProperties.hh"
 36 #include "G4NuclearLevelData.hh"
 37 #include "G4DeexPrecoParameters.hh"
 38 #include "G4VCoulombBarrier.hh"
 39 #include "G4InterfaceToXS.hh"
 40 
 41 G4VPreCompoundFragment::G4VPreCompoundFragment(
 42   const G4ParticleDefinition* part, G4VCoulombBarrier* aCoulombBarrier)
 43   : theA(part->GetBaryonNumber()),
 44     theZ(G4lrint(part->GetPDGCharge()/CLHEP::eplus)),
 45     particle(part),
 46     theCoulombBarrierPtr(aCoulombBarrier)
 47 {
 48   theMass = particle->GetPDGMass();
 49   fNucData = G4NuclearLevelData::GetInstance();
 50   theParameters = fNucData->GetParameters();
 51   OPTxs = theParameters->GetDeexModelType();
 52   g4calc = G4Pow::GetInstance();
 53 
 54   if (1 == theZ && 1 == theA) { index = 1; }
 55   else if (1 == theZ && 2 == theA) { index = 2; }
 56   else if (1 == theZ && 3 == theA) { index = 3; }
 57   else if (2 == theZ && 3 == theA) { index = 4; }
 58   else if (2 == theZ && 4 == theA) { index = 5; }
 59 
 60   if (OPTxs == 1) {
 61     fXSection = new G4InterfaceToXS(particle, index);
 62   }
 63 }
 64 
 65 G4VPreCompoundFragment::~G4VPreCompoundFragment()
 66 {
 67   delete theCoulombBarrierPtr;
 68   delete fXSection;
 69 }
 70 
 71 std::ostream& 
 72 operator << (std::ostream &out, const G4VPreCompoundFragment &theFragment)
 73 {
 74   out << &theFragment;
 75   return out; 
 76 }
 77 
 78 std::ostream& 
 79 operator << (std::ostream &out, const G4VPreCompoundFragment *theFragment)
 80 {
 81   out 
 82     << "PreCompoundModel Emitted Fragment: Z= " << theFragment->GetZ() 
 83     << " A= " << theFragment->GetA()
 84     << " Mass(GeV)= " << theFragment->GetNuclearMass()/CLHEP::GeV;
 85   return out;
 86 }
 87 
 88 G4bool 
 89 G4VPreCompoundFragment::Initialize(const G4Fragment& aFragment)
 90 {
 91   theFragA = aFragment.GetA_asInt();
 92   theFragZ = aFragment.GetZ_asInt();
 93   theResA = theFragA - theA;
 94   theResZ = theFragZ - theZ;
 95 
 96   theMinKinEnergy = theMaxKinEnergy = theCoulombBarrier = 0.0;
 97   if ((theResA < theResZ) || (theResA < theA) || (theResZ < theZ)
 98       || (theResA == theA && theResZ < theZ)
 99       || ((theResA > 1) && (theResA == theResZ || theResZ == 0))) {
100     return false;
101   }
102   theResMass = G4NucleiProperties::GetNuclearMass(theResA, theResZ);
103   G4double Ecm = aFragment.GetMomentum().m();
104   if (Ecm <= theResMass + theMass) { return 0.0; }
105 
106   theResA13 = g4calc->Z13(theResA);
107 
108   G4double elim = 0.0;
109   if (0 < theZ) {
110     theCoulombBarrier = theCoulombBarrierPtr->
111       GetCoulombBarrier(theResA, theResZ, aFragment.GetExcitationEnergy());
112     elim = (0 < OPTxs) ? theCoulombBarrier*0.5 : theCoulombBarrier;
113   }
114       
115   // Compute Maximal Kinetic Energy which can be carried by fragments 
116   // after separation - the true assimptotic value
117   theMaxKinEnergy =
118     0.5*((Ecm - theResMass)*(Ecm + theResMass) + theMass*theMass)/Ecm - theMass;
119   G4double resM = Ecm - theMass - elim;
120   if (resM < theResMass) { return false; }
121   theMinKinEnergy =
122     0.5*((Ecm - resM)*(Ecm + resM) + theMass*theMass)/Ecm - theMass;
123 
124   if (theMinKinEnergy >= theMaxKinEnergy) { return false; }
125   // Calculate masses
126   theReducedMass = theResMass*theMass/(theResMass + theMass);
127 
128   // Compute Binding Energies for fragments 
129   // needed to separate a fragment from the nucleus
130   theBindingEnergy = theResMass + theMass - aFragment.GetGroundStateMass();
131   return true;
132 }
133