Geant4-11
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//#define debugGammaParticipants
35
37{
38 // Check reaction threshold - goes to CheckThreshold
39
42
43 G4LorentzVector aPrimaryMomentum(thePrimary.GetMomentum(), thePrimary.GetTotalEnergy());
44 G4LorentzVector aTargetNMomentum(0.,0.,0.,938.);
45 if((!(aPrimaryMomentum.e()>-1)) && (!(aPrimaryMomentum.e()<1)) )
46 {
47 throw G4HadronicException(__FILE__, __LINE__,
48 "G4GammaParticipants::SelectInteractions: primary nan energy.");
49 }
50 G4double S = (aPrimaryMomentum + aTargetNMomentum).mag2();
51 G4double ThresholdMass = thePrimary.GetMass() + 938.;
53
54 #ifdef debugGammaParticipants
55 G4cout <<G4endl<< "Gamma Participants - SelectInteractions " << G4endl;
56 G4cout << "Energy and Nucleus Mass N "<<thePrimary.GetTotalEnergy()<<" "<<theNucleus->GetMassNumber()<<G4endl;
57 G4cout << "SqrtS ThresholdMass ModelMode " <<std::sqrt(S)<<" "<<ThresholdMass<<" "<<ModelMode<< G4endl;
58 G4cout << "ThresholdParameter QGSMThreshold "<<ThresholdParameter<<" "<<QGSMThreshold<<G4endl;
59 #endif
60
61
62 if (sqr(ThresholdMass + ThresholdParameter) > S)
63 {
65 }
66
67 if (sqr(ThresholdMass + QGSMThreshold) > S)
68 {
70 }
71
72 #ifdef debugGammaParticipants
73 G4cout << "Interaction type (ModelMode) 0 - SOFT, 1 - DIFFRACTIVE: "<<ModelMode<< G4endl;
74 #endif
75
76 std::for_each(theInteractions.begin(), theInteractions.end(), DeleteInteractionContent());
77 theInteractions.clear();
78
79// #ifdef debug_G4GammaParticipants
80// G4double eK = thePrimary.GetKineticEnergy()/GeV;
81// G4int nucleonCount = theNucleus->GetMassNumber();
82// #endif
83
85 G4int NucleonNo=0;
86
88 G4Nucleon * pNucleon =0;
89
90 while( (pNucleon = theNucleus->GetNextNucleon()) ) {if(NucleonNo == theCurrent) break; NucleonNo++;}
91
92 if ( pNucleon ) {
93
94 G4QGSMSplitableHadron* aTarget = new G4QGSMSplitableHadron(*pNucleon);
95 pNucleon->Hit(aTarget);
96
97 if( (0.06 > G4UniformRand() &&(ModelMode==SOFT)) || (ModelMode==DIFFRACTIVE ) ) // (false) //
98 { // Diffractive interaction
101
102 aInteraction->SetTarget(aTarget);
103 aInteraction->SetTargetNucleon(pNucleon);
104 aTarget->SetCollisionCount(0);
105 aTarget->SetStatus(1); // Mark that is Diffr. interaction
106
107 aInteraction->SetNumberOfDiffractiveCollisions(1);
108 aInteraction->SetNumberOfSoftCollisions(0);
109 aInteraction->SetStatus(1);
110
111 theInteractions.push_back(aInteraction);
112 }
113 else
114 {
115 // nondiffractive soft interaction occurs
116 aTarget->IncrementCollisionCount(1);
117 aTarget->SetStatus(0);
118 theTargets.push_back(aTarget);
119
122
123 G4InteractionContent * aInteraction =
125 aInteraction->SetTarget(aTarget);
126 aInteraction->SetTargetNucleon(pNucleon);
127 aInteraction->SetNumberOfSoftCollisions(1);
128 aInteraction->SetStatus(0); // Mark that is non-Diffr. interaction
129 theInteractions.push_back(aInteraction);
130 }
131 }
133}
G4double S(G4double temp)
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define TRUE
Definition: Globals.hh:27
#define G4UniformRand()
Definition: Randomize.hh:52
virtual G4VSplitableHadron * SelectInteractions(const G4ReactionProduct &thePrimary)
void SetNumberOfDiffractiveCollisions(int)
void SetTargetNucleon(G4Nucleon *aNucleon)
void SetTarget(G4VSplitableHadron *aTarget)
void SetStatus(G4int aValue)
void Hit(G4VSplitableHadron *aHit)
Definition: G4Nucleon.hh:91
const G4double ThresholdParameter
std::vector< G4InteractionContent * > theInteractions
G4QGSMSplitableHadron * theProjectileSplitable
const G4double QGSMThreshold
std::vector< G4VSplitableHadron * > theTargets
G4double GetTotalEnergy() const
G4ThreeVector GetMomentum() const
G4double GetMass() const
virtual G4Nucleon * GetNextNucleon()=0
virtual G4bool StartLoop()=0
virtual G4int GetMassNumber()=0
G4V3DNucleus * theNucleus
void SetStatus(const G4int aStatus)
void SetCollisionCount(G4int aCount)
void IncrementCollisionCount(G4int aCount)
T sqr(const T &x)
Definition: templates.hh:128