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