Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4GammaParticipants.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#include "globals.hh"
28#include "G4LorentzVector.hh"
29#include "G4V3DNucleus.hh"
30#include <utility>
31
32// Class G4GammaParticipants
33
34// J.P. Wellisch, April 2002
35// new participants class for gamma nuclear, with this design more can come with
36// cross-section based, and quasi-eiconal model based modelling
37//
38// 20110805 M. Kelsey -- Follow change to G4V3DNucleus::GetNucleons()
39
40G4VSplitableHadron* G4GammaParticipants::SelectInteractions(const G4ReactionProduct &thePrimary)
41{
42 // Check reaction threshold - goes to CheckThreshold
43 G4VSplitableHadron* aProjectile = new G4QGSMSplitableHadron(thePrimary, TRUE); // @@@ check the TRUE
44
45 const std::vector<G4Nucleon>& theTargetNuc = theNucleus->GetNucleons();
46 G4LorentzVector aPrimaryMomentum(thePrimary.GetMomentum(), thePrimary.GetTotalEnergy());
47 if((!(aPrimaryMomentum.e()>-1)) && (!(aPrimaryMomentum.e()<1)) )
48 {
49 throw G4HadronicException(__FILE__, __LINE__,
50 "G4GammaParticipants::SelectInteractions: primary nan energy.");
51 }
52 G4double S = (aPrimaryMomentum + theTargetNuc[0].Get4Momentum()).mag2();
53 G4double ThresholdMass = thePrimary.GetMass() + theTargetNuc[0].GetDefinition()->GetPDGMass();
55 if (sqr(ThresholdMass + ThresholdParameter) > S)
56 {
58 //throw G4HadronicException(__FILE__, __LINE__, "Initial energy is too low. The 4-vectors of the input are inconsistant with the particle masses.");
59 }
60 if (sqr(ThresholdMass + QGSMThreshold) > S) // thus only diffractive in cascade!
61 {
63 }
64
65 // first find the collisions HPW
66 std::for_each(theInteractions.begin(), theInteractions.end(), DeleteInteractionContent());
67 theInteractions.clear();
68 G4int totalCuts = 0;
69
70 #ifdef debug_G4GammaParticipants
71 G4double eK = thePrimary.GetKineticEnergy()/GeV;
72 G4int nucleonCount = theTargetNuc.size(); // debug
73 #endif
74
75 G4int theCurrent = static_cast<G4int> (theTargetNuc.size()*G4UniformRand());
76 const G4Nucleon& pNucleon = theTargetNuc[theCurrent];
77 G4QGSMSplitableHadron* aTarget = new G4QGSMSplitableHadron(pNucleon);
78 theTargets.push_back(aTarget);
79 const_cast<G4Nucleon&>(pNucleon).Hit(aTarget);
80 if ( (0.06 > G4UniformRand() &&(ModelMode==SOFT)) || (ModelMode==DIFFRACTIVE ) )
81 {
82 // diffractive interaction occurs
84 {
85 theSingleDiffExcitation.ExciteParticipants(aProjectile, aTarget);
86 }
87 else
88 {
89 theDiffExcitaton.ExciteParticipants(aProjectile, aTarget);
90 }
91 G4InteractionContent * aInteraction = new G4InteractionContent(aProjectile);
92 aInteraction->SetTarget(aTarget);
93 theInteractions.push_back(aInteraction);
94 aInteraction->SetNumberOfDiffractiveCollisions(1);
95 totalCuts += 1;
96 }
97 else
98 {
99 // nondiffractive soft interaction occurs
100 aTarget->IncrementCollisionCount(1);
101 aProjectile->IncrementCollisionCount(1);
102 G4InteractionContent * aInteraction = new G4InteractionContent(aProjectile);
103 aInteraction->SetTarget(aTarget);
104 aInteraction->SetNumberOfSoftCollisions(1);
105 theInteractions.push_back(aInteraction);
106 totalCuts += 1;
107 }
108 return aProjectile;
109}
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4UniformRand()
Definition: Randomize.hh:53
void SetNumberOfDiffractiveCollisions(int)
void SetTarget(G4VSplitableHadron *aTarget)
virtual G4bool ExciteParticipants(G4VSplitableHadron *aPartner, G4VSplitableHadron *bPartner) const
const G4double ThresholdParameter
G4SingleDiffractiveExcitation theSingleDiffExcitation
std::vector< G4InteractionContent * > theInteractions
G4QGSDiffractiveExcitation theDiffExcitaton
const G4double QGSMThreshold
std::vector< G4VSplitableHadron * > theTargets
G4double GetKineticEnergy() const
G4double GetTotalEnergy() const
G4ThreeVector GetMomentum() const
G4double GetMass() const
G4bool ExciteParticipants(G4VSplitableHadron *aPartner, G4VSplitableHadron *bPartner) const
virtual const std::vector< G4Nucleon > & GetNucleons()=0
G4V3DNucleus * theNucleus
void IncrementCollisionCount(G4int aCount)
#define TRUE
Definition: globals.hh:55
T sqr(const T &x)
Definition: templates.hh:145