Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4InuclEvaporation.cc
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: G4InuclEvaporation.cc 66241 2012-12-13 18:34:42Z gunter $
27 //
28 // 20100114 M. Kelsey -- Remove G4CascadeMomentum, use G4LorentzVector directly
29 // 20100405 M. Kelsey -- Pass const-ref std::vector<>
30 // 20100413 M. Kelsey -- Pass G4CollisionOutput by ref to ::collide(), use
31 // const_iterator.
32 // 20100428 M. Kelsey -- Use G4InuclParticleNames enum
33 // 20100429 M. Kelsey -- Change "case gamma:" to "case photon:"
34 // 20100517 M. Kelsey -- Follow new ctors for G4*Collider family. Make
35 // G4EvaporationInuclCollider a data member.
36 // 20100520 M. Kelsey -- Simplify collision loop, move momentum rotations to
37 // G4CollisionOutput, copy G4DynamicParticle directly from
38 // G4InuclParticle, no switch-block required. Fix scaling factors.
39 // 20100914 M. Kelsey -- Migrate to integer A and Z
40 // 20100924 M. Kelsey -- Migrate to "OutgoingNuclei" names in CollisionOutput
41 // 20110728 M. Kelsey -- Fix Coverity #28776, remove return after throw.
42 // 20120608 M. Kelsey -- Fix variable-name "shadowing" compiler warnings.
43 // 20140310 M. Kelsey -- *TEMPORARY* const-cast G4PD* for G4Fragment ctor.
44 
45 #include <numeric>
46 
47 #include "G4InuclEvaporation.hh"
48 #include "globals.hh"
49 #include "G4SystemOfUnits.hh"
50 #include "G4IonTable.hh"
51 #include "G4V3DNucleus.hh"
54 #include "G4InuclNuclei.hh"
55 #include "G4Track.hh"
56 #include "G4Nucleus.hh"
57 #include "G4Nucleon.hh"
58 #include "G4NucleiModel.hh"
59 #include "G4HadronicException.hh"
60 #include "G4LorentzVector.hh"
63 #include "G4InuclParticle.hh"
64 #include "G4CollisionOutput.hh"
65 #include "G4InuclParticleNames.hh"
66 
67 using namespace G4InuclParticleNames;
68 
69 typedef std::vector<G4InuclElementaryParticle>::const_iterator particleIterator;
70 typedef std::vector<G4InuclNuclei>::const_iterator nucleiIterator;
71 
72 
74  : verboseLevel(0), evaporator(new G4EvaporationInuclCollider) {}
75 
77  throw G4HadronicException(__FILE__, __LINE__, "G4InuclEvaporation::copy_constructor meant to not be accessable.");
78 }
79 
81  delete evaporator;
82 }
83 
84 const G4InuclEvaporation & G4InuclEvaporation::operator=(const G4InuclEvaporation &) {
85  throw G4HadronicException(__FILE__, __LINE__, "G4InuclEvaporation::operator= meant to not be accessable.");
86 }
87 
88 
89 G4bool G4InuclEvaporation::operator==(const G4InuclEvaporation &) const {
90  return false;
91 }
92 
93 G4bool G4InuclEvaporation::operator!=(const G4InuclEvaporation &) const {
94  return true;
95 }
96 
98  verboseLevel = verbose;
99 }
100 
102  G4FragmentVector* theResult = new G4FragmentVector;
103 
104  if (theNucleus.GetExcitationEnergy() <= 0.0) { // Check that Excitation Energy > 0
105  theResult->push_back(new G4Fragment(theNucleus));
106  return theResult;
107  }
108 
109  G4int A = theNucleus.GetA_asInt();
110  G4int Z = theNucleus.GetZ_asInt();
111  G4double mTar = G4NucleiProperties::GetNuclearMass(A, Z); // Mass of the target nucleus
112 
113  G4ThreeVector momentum = theNucleus.GetMomentum().vect();
114  G4double exitationE = theNucleus.GetExcitationEnergy();
115 
116  G4double mass = mTar;
117  G4ThreeVector boostToLab( momentum/mass );
118 
119  if ( verboseLevel > 2 )
120  G4cout << " G4InuclEvaporation : initial kinematics : boostToLab vector = " << boostToLab << G4endl
121  << " excitation energy : " << exitationE << G4endl;
122 
123  if (verboseLevel > 2) {
124  G4cout << "G4InuclEvaporation::BreakItUp >>> A: " << A << " Z: " << Z
125  << " exitation E: " << exitationE << " mass: " << mTar/GeV << " GeV"
126  << G4endl;
127  };
128 
129  G4InuclNuclei* nucleus = new G4InuclNuclei(A, Z);
130  nucleus->setExitationEnergy(exitationE);
131 
132  G4CollisionOutput output;
133  evaporator->collide(0, nucleus, output);
134 
135  const std::vector<G4InuclNuclei>& outgoingNuclei = output.getOutgoingNuclei();
136  const std::vector<G4InuclElementaryParticle>& particles = output.getOutgoingParticles();
137 
138  G4double eTot=0.0;
139  G4int i=1;
140 
141  if (!particles.empty()) {
142  G4int outgoingType;
143  particleIterator ipart = particles.begin();
144  for (; ipart != particles.end(); ipart++) {
145  outgoingType = ipart->type();
146 
147  if (verboseLevel > 2) {
148  G4cout << "Evaporated particle: " << i << " of type: "
149  << outgoingType << G4endl;
150  i++;
151  }
152 
153  eTot += ipart->getEnergy();
154 
155  G4LorentzVector vlab = ipart->getMomentum().boost(boostToLab);
156 
157  // TEMPORARY: Remove constness on PD until G4Fragment is fixed
158  theResult->push_back( new G4Fragment(vlab, ipart->getDefinition()) );
159  }
160  }
161 
162  if (!outgoingNuclei.empty()) {
163  nucleiIterator ifrag = outgoingNuclei.begin();
164  for (i=1; ifrag != outgoingNuclei.end(); ifrag++) {
165  if (verboseLevel > 2) {
166  G4cout << " Nuclei fragment: " << i << G4endl; i++;
167  }
168 
169  eTot += ifrag->getEnergy();
170 
171  G4LorentzVector vlab = ifrag->getMomentum().boost(boostToLab);
172 
173  G4int fragA = ifrag->getA();
174  G4int fragZ = ifrag->getZ();
175  if (verboseLevel > 2) {
176  G4cout << "boosted v" << vlab << G4endl;
177  }
178  theResult->push_back( new G4Fragment(fragA, fragZ, vlab) );
179  }
180  }
181 
182  return theResult;
183 }
static G4double GetNuclearMass(const G4double A, const G4double Z)
void setVerboseLevel(const G4int verbose)
G4FragmentVector * BreakItUp(const G4Fragment &theNucleus)
int G4int
Definition: G4Types.hh:78
void setExitationEnergy(G4double e)
Hep3Vector vect() const
G4GLOB_DLL std::ostream G4cout
double A(double temperature)
G4int GetA_asInt() const
Definition: G4Fragment.hh:266
bool G4bool
Definition: G4Types.hh:79
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:307
HepLorentzVector & boost(double, double, double)
std::vector< G4Fragment * > G4FragmentVector
Definition: G4Fragment.hh:63
virtual void collide(G4InuclParticle *bullet, G4InuclParticle *target, G4CollisionOutput &globalOutput)
const std::vector< G4InuclNuclei > & getOutgoingNuclei() const
const std::vector< G4InuclElementaryParticle > & getOutgoingParticles() const
static constexpr double GeV
Definition: G4SIunits.hh:217
G4int GetZ_asInt() const
Definition: G4Fragment.hh:271
std::vector< G4InuclNuclei >::const_iterator nucleiIterator
std::vector< G4InuclElementaryParticle >::iterator particleIterator
Definition: G4BigBanger.cc:65
#define G4endl
Definition: G4ios.hh:61
double getZ() const
double G4double
Definition: G4Types.hh:76
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:283