Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4LightTargetCollider.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// //
28// File: G4LightTargetCollider.cc //
29// Date: 30 September 2019 //
30// Author: Dennis Wright (SLAC) //
31// //
32// Description: model for collision of elementary particles with light //
33// targets (H, D, T, 3He) //
34// //
35////////////////////////////////////////////////////////////////////////////////
36
38#include "G4CascadeChannel.hh"
41#include "G4CollisionOutput.hh"
44#include "G4InuclNuclei.hh"
45#include "G4NucleiModel.hh"
46#include "G4LorentzConvertor.hh"
47#include "G4Deuteron.hh"
48#include "G4Gamma.hh"
49#include "G4PionZero.hh"
50#include "G4PionPlus.hh"
51#include "G4PionMinus.hh"
52#include "G4RandomDirection.hh"
53
54
56 : G4CascadeColliderBase("G4LightTargetCollider"),
57 theElementaryParticleCollider(new G4ElementaryParticleCollider)
58{
59 mP = G4Proton::Proton()->GetPDGMass()/CLHEP::GeV;
60 mN = G4Neutron::Neutron()->GetPDGMass()/CLHEP::GeV;
61 mD = G4Deuteron::Deuteron()->GetPDGMass()/CLHEP::GeV;
62 pFermiD = 0.045; // Fermi momentum of nucleon in deuteron Hulthen potential
63}
64
66 delete theElementaryParticleCollider;
67}
68
69
70// Set verbosity and pass on to member objects
76
77
79 G4InuclParticle* target,
80 G4CollisionOutput& globalOutput)
81{
82 if (verboseLevel) {
83 G4cout << " >>> G4LightTargetCollider::collide" << G4endl;
84 G4cout << " Projectile: " << bullet->getDefinition()->GetParticleName() << G4endl;
85 G4cout << " Target: " << target->getDefinition()->GetParticleName() << G4endl;
86 }
87
88 G4double ke = bullet->getKineticEnergy();
89
90 if (target->getDefinition() == G4Proton::Proton() ) {
91 if (ke < 0.1447) {
92 // Below threshold lab energy for pi0 creation
93 globalOutput.trivialise(bullet, target);
94 } else {
95 // Need inelastic cross section in this class if we want to replace ElementaryParticleCollider
96 // with SingleNucleonScattering
97 theElementaryParticleCollider->collide(bullet, target, globalOutput);
98 if (globalOutput.numberOfOutgoingParticles() == 0) globalOutput.trivialise(bullet, target);
99 }
100
101 } else if (target->getDefinition() == G4Deuteron::Deuteron() ) {
102
103 if (ke < mP + mN - mD) {
104 // Should not happen as long as inelastic cross section is zero
105 G4Exception("G4LightTargetCollider::collide()","HAD_BERT_201",
106 JustWarning, "Projectile energy below reaction threshold");
107 globalOutput.trivialise(bullet, target);
108
109 } else {
110 // Get p, n and deuteron cross sections; use lab energy to access
113 G4double gammaDXS = GammaDCrossSection(ke);
114
115 G4double probP = 0.0;
116 G4double probN = 0.0;
117 // Highest threshold is 0.152 (for gamma p -> n pi+)
118 // Because of Fermi momentum in deuteron, raise this to 0.159
119 if (ke > 0.159) {
120 G4double totalDXS = gammaPXS + gammaNXS + gammaDXS;
121 probP = gammaPXS/totalDXS;
122 probN = (gammaPXS+gammaNXS)/totalDXS;
123 }
124
125 G4double rndm = G4UniformRand();
126 if (rndm < probP) {
127 // Generate Fermi momenta of bullet and target
128 G4ThreeVector fermiMomentum = pFermiD*G4RandomDirection();
129 G4LorentzVector protonMomentum(fermiMomentum, std::sqrt(mP*mP + pFermiD*pFermiD) );
130 G4LorentzVector neutronMomentum(-fermiMomentum, std::sqrt(mN*mN + pFermiD*pFermiD) );
131
132 G4LorentzVector bulletMomentum = bullet->getMomentum();
133 G4ThreeVector betacm = bulletMomentum.findBoostToCM(protonMomentum);
134
135 // First boost bullet and target so that target is at rest
136 G4ThreeVector toProtonRest = -protonMomentum.boostVector();
137 protonMomentum.boost(toProtonRest);
138 bulletMomentum.boost(toProtonRest);
139
140 G4InuclElementaryParticle projectile(bulletMomentum, bullet->getDefinition() );
141 G4InuclElementaryParticle targetNucleon(protonMomentum, G4Proton::Proton() );
142 G4InuclElementaryParticle spectatorNucleon(neutronMomentum, G4Neutron::Neutron() );
143 ScatteringProducts products = SingleNucleonScattering(projectile, targetNucleon);
144
145 // Particles from SingleNucleonScattering are in CM frame of projectile
146 // and moving proton. Transform back to lab frame with -betacm, then
147 // add them to outgoing list.
148 globalOutput.reset();
149 G4LorentzVector temp;
150 for (G4int i = 0; i < G4int(products.size()); i++) {
151 temp = products[i].getMomentum();
152 temp.boost(-betacm);
153 products[i].setMomentum(temp);
154 globalOutput.addOutgoingParticle(products[i]);
155 }
156
157 // Add the recoil nucleon unmodified
158 globalOutput.addOutgoingParticle(spectatorNucleon);
159
160 } else if (rndm < probN) {
161 G4ThreeVector fermiMomentum = pFermiD*G4RandomDirection();
162 G4LorentzVector protonMomentum(fermiMomentum, std::sqrt(mP*mP + pFermiD*pFermiD) );
163 G4LorentzVector neutronMomentum(-fermiMomentum, std::sqrt(mN*mN + pFermiD*pFermiD) );
164
165 G4LorentzVector bulletMomentum = bullet->getMomentum();
166 G4ThreeVector betacm = bulletMomentum.findBoostToCM(neutronMomentum);
167
168 // First boost bullet and target so that target is at rest
169 G4ThreeVector toNeutronRest = -neutronMomentum.boostVector();
170 neutronMomentum.boost(toNeutronRest);
171 bulletMomentum.boost(toNeutronRest);
172
173 G4InuclElementaryParticle projectile(bulletMomentum, bullet->getDefinition() );
174 G4InuclElementaryParticle targetNucleon(neutronMomentum, G4Neutron::Neutron() );
175 G4InuclElementaryParticle spectatorNucleon(protonMomentum, G4Proton::Proton() );
176
177 ScatteringProducts products = SingleNucleonScattering(projectile, targetNucleon);
178
179 // Particles from SingleNucleonScattering are in CM frame of projectile
180 // and moving neutron. Transform back to lab frame with -betacm, then add
181 // them to outgoing list
182 globalOutput.reset();
183 G4LorentzVector temp;
184 for (G4int i = 0; i < G4int(products.size()); i++) {
185 temp = products[i].getMomentum();
186 temp.boost(-betacm);
187 products[i].setMomentum(temp);
188 globalOutput.addOutgoingParticle(products[i]);
189 }
190
191 // Add the recoil nucleon unmodified
192 globalOutput.addOutgoingParticle(spectatorNucleon);
193
194 } else {
195 NucleonPair products = AbsorptionOnDeuteron(bullet);
196 globalOutput.reset();
197 globalOutput.addOutgoingParticle(products.first);
198 globalOutput.addOutgoingParticle(products.second);
199 }
200 } // Energy above threshold ?
201
202 // Test code
203 // G4int numPart = globalOutput.numberOfOutgoingParticles();
204 // std::vector<G4InuclElementaryParticle> testList = globalOutput.getOutgoingParticles();
205 // G4LorentzVector sumP;
206 // G4cout << " Global output " << G4endl;
207 // for (G4int i = 0; i < numPart; i++) {
208 // sumP += testList[i].getMomentum();
209 // G4cout << testList[i] << G4endl;
210 // }
211 // G4cout << " Global 4-momentum sum = " << sumP << G4endl;
212 // G4cout << " Initial lab energy = " << mD + bullet->getEnergy() << G4endl;
213
214 } else {
215 G4Exception("G4LightTargetCollider::collide()","HAD_BERT_203",
216 FatalException, "Scattering from this target not implemented");
217 }
218
219 return;
220}
221
222
223G4double G4LightTargetCollider::GammaDCrossSection(G4double gammaEnergy)
224{
225 // Gamma deuteron cross section in mb parameterized from JLab data
226 // No parameterization needed below pi0 threshold where cross section
227 // is 100% disintegration
228 G4double sigma = 1000.0;
229 G4double term = 0.;
230 if (gammaEnergy > 0.144 && gammaEnergy < 0.42) {
231 term = (gammaEnergy - 0.24)/0.155;
232 sigma = 0.065*std::exp(-term*term);
233 } else if (gammaEnergy >= 0.42) {
234 sigma = 0.000526/gammaEnergy/gammaEnergy/gammaEnergy/gammaEnergy;
235 }
236
237 return sigma;
238}
239
240
241NucleonPair G4LightTargetCollider::AbsorptionOnDeuteron(G4InuclParticle* bullet)
242{
243 // Do break-up in center of mass, convert to lab frame before returning
244 // particles
245
246 G4double bulletMass = bullet->getMass();
247 G4double bulletE = bullet->getEnergy();
248
249 G4double S = bulletMass*bulletMass + mD*mD + 2.*mD*bulletE;
250
251 G4double qcm = 0.;
252 G4int outType1 = 0;
253 G4int outType2 = 0;
254 G4LorentzVector Mom1;
255 G4LorentzVector Mom2;
256
257 // Set up outgoing particle types
258 if (bullet->getDefinition() == G4Gamma::Gamma() ||
259 bullet->getDefinition() == G4PionZero::PionZero() ) {
260 qcm = std::sqrt( (S - (mP + mN)*(mP + mN)) * (S - (mP - mN)*(mP - mN))/S/4.);
261 Mom1.setE(std::sqrt(mP*mP + qcm*qcm) );
263 Mom2.setE(std::sqrt(mN*mN + qcm*qcm) );
265
266 } else if (bullet->getDefinition() == G4PionPlus::PionPlus() ) {
267 qcm = std::sqrt( (S - 4.*mP*mP)/4.);
268 Mom1.setE(std::sqrt(mP*mP + qcm*qcm) );
270 Mom2.setE(std::sqrt(mP*mP + qcm*qcm) );
272
273 } else if (bullet->getDefinition() == G4PionMinus::PionMinus() ) {
274 qcm = std::sqrt( (S - 4.*mN*mN)/4.);
275 Mom1.setE(std::sqrt(mN*mN + qcm*qcm) );
277 Mom2.setE(std::sqrt(mN*mN + qcm*qcm) );
279
280 } else {
281 G4Exception("G4LightTargetCollider::collide()","HAD_BERT_204",
282 FatalException, "Illegal bullet type");
283 }
284
285 // Sample angular distribution, assuming 100% S wave (no D-wave)
286 G4ThreeVector qVect = qcm*G4RandomDirection();
287 Mom1.setVect(qVect);
288 Mom2.setVect(-qVect);
289
290 // Boost to lab frame
291 G4ThreeVector betacm(0., 0., bullet->getMomModule()/(bulletE + mD) );
292 Mom1.boost(betacm);
293 Mom2.boost(betacm);
294
295 G4InuclElementaryParticle particle1(Mom1, outType1);
296 G4InuclElementaryParticle particle2(Mom2, outType2);
297 NucleonPair nucleon_pair(particle1, particle2);
298
299 // if pion, use parameterization of B.G. Ritchie, PRC 44, 533 (1991)
300 // Total cross section: 1/E + Lorentzian
301
302 return nucleon_pair;
303}
304
305
307G4LightTargetCollider::SingleNucleonScattering(const G4InuclElementaryParticle& projectile,
308 const G4InuclElementaryParticle& nucleon)
309{
310 // At this point projectile and nucleon momenta are in nucleon rest frame
311 G4int reactionIndex = G4InuclElementaryParticle::type(projectile.getDefinition() )
312 * G4InuclElementaryParticle::type(nucleon.getDefinition() );
313
314 const G4CascadeChannel* xsecTable = G4CascadeChannelTables::GetTable(reactionIndex);
315 G4double ke = projectile.getKineticEnergy();
316 G4int mult = xsecTable->getMultiplicity(ke);
317
318 std::vector<G4double> masses;
319 G4double mass = 0.0;
320 G4LorentzVector totalMom = projectile.getMomentum() + nucleon.getMomentum();
321 G4double Ecm = totalMom.mag();
322
323 std::vector<G4LorentzVector> cmMomenta;
324 std::vector<G4int> particle_kinds;
325 G4int itry = 0;
326 G4int itry_max = 200;
327 G4bool generate = true;
328
329 while (mult > 1) {
330 itry = 0;
331 generate = true;
332 while (generate && itry < itry_max) {
333 particle_kinds.clear();
334 xsecTable->getOutgoingParticleTypes(particle_kinds, mult, ke);
335 masses.clear();
336 for (G4int i = 0; i < mult; i++) {
337 mass = G4InuclElementaryParticle::getParticleMass(particle_kinds[i]);
338 masses.push_back(mass);
339 }
340
341 fsGen.Configure(const_cast<G4InuclElementaryParticle*>(&projectile),
342 const_cast<G4InuclElementaryParticle*>(&nucleon),
343 particle_kinds);
344 // Generate final state in CM of projectile and at-rest nucleon
345 cmMomenta.clear();
346 generate = !fsGen.Generate(Ecm, masses, cmMomenta);
347 itry++;
348 } // while
349
350 if (itry == itry_max) mult--;
351 else break;
352
353 } // while mult
354
355 ScatteringProducts finalState;
356 if (mult < 2) {
357 G4Exception("G4LightTargetCollider::SingleNucleonScattering()","HAD_BERT_202",
358 JustWarning, "Failed to generate final state");
359 // Final state particles not in CM - just using them as dummies
360 finalState.push_back(projectile);
361 finalState.push_back(nucleon);
362
363 } else {
364 for (G4int i = 0; i < mult; i++) {
365 G4InuclElementaryParticle fsPart(cmMomenta[i], particle_kinds[i]);
366 finalState.push_back(fsPart);
367 }
368 }
369
370 return finalState;
371}
G4double S(G4double temp)
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::vector< G4InuclElementaryParticle > ScatteringProducts
std::pair< G4InuclElementaryParticle, G4InuclElementaryParticle > NucleonPair
G4ThreeVector G4RandomDirection()
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition Randomize.hh:52
Hep3Vector boostVector() const
HepLorentzVector & boost(double, double, double)
void setVect(const Hep3Vector &)
Hep3Vector findBoostToCM() const
static const G4CascadeChannel * GetTable(G4int initialState)
virtual G4int getMultiplicity(G4double ke) const =0
virtual G4double getCrossSection(double ke) const =0
virtual void getOutgoingParticleTypes(std::vector< G4int > &kinds, G4int mult, G4double ke) const =0
virtual void setVerboseLevel(G4int verbose=0)
void Configure(G4InuclElementaryParticle *bullet, G4InuclElementaryParticle *target, const std::vector< G4int > &particle_kinds)
G4int numberOfOutgoingParticles() const
void addOutgoingParticle(const G4InuclElementaryParticle &particle)
void setVerboseLevel(G4int verbose)
void trivialise(G4InuclParticle *bullet, G4InuclParticle *target)
static G4Deuteron * Deuteron()
Definition G4Deuteron.cc:90
void collide(G4InuclParticle *bullet, G4InuclParticle *target, G4CollisionOutput &output)
static G4Gamma * Gamma()
Definition G4Gamma.cc:81
G4bool Generate(G4double initialMass, const std::vector< G4double > &masses, std::vector< G4LorentzVector > &finalState)
static G4double getParticleMass(G4int type)
const G4ParticleDefinition * getDefinition() const
G4double getKineticEnergy() const
G4double getMass() const
G4LorentzVector getMomentum() const
G4double getMomModule() const
G4double getEnergy() const
void setVerboseLevel(G4int verbose=0)
void collide(G4InuclParticle *bullet, G4InuclParticle *target, G4CollisionOutput &globalOutput)
static G4Neutron * Neutron()
Definition G4Neutron.cc:101
const G4String & GetParticleName() const
static G4PionMinus * PionMinus()
static G4PionPlus * PionPlus()
Definition G4PionPlus.cc:93
static G4PionZero * PionZero()
static G4Proton * Proton()
Definition G4Proton.cc:90
void generate(const G4double sqrtS, ParticleList &particles)
Generate an event in the CM system.
G4bool nucleon(G4int ityp)