source: trunk/source/processes/hadronic/models/pre_equilibrium/exciton_model/src/G4VPreCompoundFragment.cc @ 829

Last change on this file since 829 was 819, checked in by garnier, 16 years ago

import all except CVS

File size: 6.2 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// $Id: G4VPreCompoundFragment.cc,v 1.4 2007/07/23 09:56:40 ahoward Exp $
28// GEANT4 tag $Name:  $
29//
30// by V. Lara
31 
32#include "G4VPreCompoundFragment.hh"
33#include "G4PreCompoundParameters.hh"
34
35G4VPreCompoundFragment::
36G4VPreCompoundFragment(const G4VPreCompoundFragment & right)
37{
38  theA = right.theA;
39  theZ = right.theZ;
40  theRestNucleusA = right.theRestNucleusA;
41  theRestNucleusZ = right.theRestNucleusZ;
42  theCoulombBarrier = right.theCoulombBarrier;
43  theCoulombBarrierPtr = right.theCoulombBarrierPtr;
44  theMaximalKineticEnergy = right.theMaximalKineticEnergy;
45  theEmissionProbability = right.theEmissionProbability;
46  theMomentum = right.theMomentum;
47  theFragmentName = right.theFragmentName;
48  theStage = right.theStage;
49}
50
51G4VPreCompoundFragment::
52G4VPreCompoundFragment(const G4double anA,
53                       const G4double aZ, 
54                       G4VCoulombBarrier* aCoulombBarrier,
55                       const G4String & aName):
56  theA(anA),theZ(aZ), 
57  theRestNucleusA(0.0),theRestNucleusZ(0.0),theCoulombBarrier(0.0),
58  theCoulombBarrierPtr(aCoulombBarrier),
59  theBindingEnergy(0.0), theMaximalKineticEnergy(-1.0),
60  theEmissionProbability(0.0), theMomentum(0.0,0.0,0.0,0.0),
61  theFragmentName(aName),theStage(0)
62{}
63
64
65
66G4VPreCompoundFragment::~G4VPreCompoundFragment()
67{
68}
69
70
71const G4VPreCompoundFragment & G4VPreCompoundFragment::
72operator= (const G4VPreCompoundFragment & right)
73{
74  if (this != &right) {
75    theA = right.theA;
76    theZ = right.theZ;
77    theRestNucleusA = right.theRestNucleusA;
78    theRestNucleusZ = right.theRestNucleusZ;
79    theCoulombBarrier = right.theCoulombBarrier;
80    theCoulombBarrierPtr = right.theCoulombBarrierPtr;
81    theMaximalKineticEnergy = right.theMaximalKineticEnergy;
82    theEmissionProbability = right.theEmissionProbability;
83    theMomentum = right.theMomentum;
84    theFragmentName = right.theFragmentName;
85    theStage = right.theStage;
86  }
87  return *this;
88}
89
90G4int G4VPreCompoundFragment::operator==(const G4VPreCompoundFragment & right) const
91{
92  return (this == (G4VPreCompoundFragment *) &right);
93}
94
95G4int G4VPreCompoundFragment::operator!=(const G4VPreCompoundFragment & right) const
96{
97  return (this != (G4VPreCompoundFragment *) &right);
98}
99
100
101std::ostream& 
102operator << (std::ostream &out, const G4VPreCompoundFragment &theFragment)
103{
104  out << &theFragment;
105  return out; 
106}
107
108
109std::ostream& 
110operator << (std::ostream &out, const G4VPreCompoundFragment *theFragment)
111{
112  std::ios::fmtflags old_floatfield = out.flags();
113  out.setf(std::ios::floatfield);
114   
115  out
116    << "PreCompound Model Emitted Fragment: A = " 
117    << std::setprecision(3) << theFragment->theA
118    << ", Z = " << std::setprecision(3) << theFragment->theZ;
119    out.setf(std::ios::scientific,std::ios::floatfield);
120    //   out
121    //     << ", U = " << theFragment->theExcitationEnergy/MeV
122    //     << " MeV" << endl
123    //     << "          P = ("
124    //     << theFragment->theMomentum.x()/MeV << ","
125    //     << theFragment->theMomentum.y()/MeV << ","
126    //     << theFragment->theMomentum.z()/MeV
127    //     << ") MeV   E = "
128    //     << theFragment->theMomentum.t()/MeV << " MeV";
129   
130    out.setf(old_floatfield,std::ios::floatfield);
131   
132    return out;
133}
134
135
136void G4VPreCompoundFragment::
137Initialize(const G4Fragment & aFragment)
138{
139  theRestNucleusA = aFragment.GetA() - theA;
140  theRestNucleusZ = aFragment.GetZ() - theZ;
141
142  if ((theRestNucleusA < theRestNucleusZ) ||
143      (theRestNucleusA < theA) ||
144      (theRestNucleusZ < theZ) ||
145      (aFragment.GetNumberOfCharged() < theZ)) // AH last argument from JMQ
146    {
147      // In order to be sure that emission probability will be 0.
148      theMaximalKineticEnergy = 0.0;
149      return;
150    }
151 
152 
153  // Calculate Coulomb barrier
154  theCoulombBarrier = theCoulombBarrierPtr->
155    GetCoulombBarrier(static_cast<G4int>(theRestNucleusA),static_cast<G4int>(theRestNucleusZ),
156                      aFragment.GetExcitationEnergy());
157 
158  // Compute Binding Energies for fragments
159  // (needed to separate a fragment from the nucleus)
160 
161  theBindingEnergy = G4NucleiProperties::GetMassExcess(static_cast<G4int>(theA),static_cast<G4int>(theZ)) +
162    G4NucleiProperties::GetMassExcess(static_cast<G4int>(theRestNucleusA),static_cast<G4int>(theRestNucleusZ)) -
163    G4NucleiProperties::GetMassExcess(static_cast<G4int>(aFragment.GetA()),static_cast<G4int>(aFragment.GetZ()));
164 
165  // Compute Maximal Kinetic Energy which can be carried by fragments after separation
166  G4double m = aFragment.GetMomentum().m();
167  G4double rm = GetRestNuclearMass();
168  G4double em = GetNuclearMass();
169  theMaximalKineticEnergy = ((m - rm)*(m + rm) + em*em)/(2.0*m) - em;
170  // - theCoulombBarrier;
171 
172  return;
173}
174
175
176
177
Note: See TracBrowser for help on using the repository browser.