Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ProjectileDiffractiveChannel Class Reference

#include <G4ProjectileDiffractiveChannel.hh>

Public Member Functions

 G4ProjectileDiffractiveChannel ()
 
 ~G4ProjectileDiffractiveChannel ()
 
G4double GetFraction (G4Nucleus &theNucleus, const G4DynamicParticle &thePrimary)
 
G4KineticTrackVectorScatter (G4Nucleus &theNucleus, const G4DynamicParticle &thePrimary)
 

Detailed Description

Definition at line 47 of file G4ProjectileDiffractiveChannel.hh.

Constructor & Destructor Documentation

◆ G4ProjectileDiffractiveChannel()

G4ProjectileDiffractiveChannel::G4ProjectileDiffractiveChannel ( )

Definition at line 49 of file G4ProjectileDiffractiveChannel.cc.

50{
51 theQDiffraction=G4QDiffractionRatio::GetPointer();
52}
static G4QDiffractionRatio * GetPointer()

◆ ~G4ProjectileDiffractiveChannel()

G4ProjectileDiffractiveChannel::~G4ProjectileDiffractiveChannel ( )

Definition at line 54 of file G4ProjectileDiffractiveChannel.cc.

55{}

Member Function Documentation

◆ GetFraction()

G4double G4ProjectileDiffractiveChannel::GetFraction ( G4Nucleus theNucleus,
const G4DynamicParticle thePrimary 
)

Definition at line 57 of file G4ProjectileDiffractiveChannel.cc.

59{
60 G4double ratio;
61 ratio=theQDiffraction->GetRatio(thePrimary.GetTotalMomentum(),
62 thePrimary.GetDefinition()->GetPDGEncoding(),
63 theNucleus.GetZ_asInt(), theNucleus.GetN_asInt() );
64 #ifdef debug_getFraction
65 G4cout << "G4ProjectileDiffractiveChannel::ratio " << ratio << G4endl;
66 #endif
67
68 return ratio;
69}
double G4double
Definition: G4Types.hh:64
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
G4ParticleDefinition * GetDefinition() const
G4double GetTotalMomentum() const
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
G4int GetN_asInt() const
Definition: G4Nucleus.hh:112
G4double GetRatio(G4double pIU, G4int prPDG, G4int tgZ, G4int tgN)

Referenced by G4TheoFSGenerator::ApplyYourself().

◆ Scatter()

G4KineticTrackVector * G4ProjectileDiffractiveChannel::Scatter ( G4Nucleus theNucleus,
const G4DynamicParticle thePrimary 
)

Definition at line 71 of file G4ProjectileDiffractiveChannel.cc.

73{
74
75
76 G4int A=theNucleus.GetA_asInt();
77 G4int Z=theNucleus.GetZ_asInt();
78
79 G4ParticleDefinition * pDef=thePrimary.GetDefinition();
80 #ifdef debug_scatter
81 G4cout << " ProjectileDiffractive: A, Z, proj-pdg" <<" "<< A <<" "<<Z
82 << " "<<" " << pDef->GetParticleName()<< G4endl;
83 #endif
84
85 G4QHadronVector * result(0);
86 G4KineticTrackVector * ktv(0);
87 G4bool tryAgain;
88
89 do {
90 tryAgain = false;
91 result=theQDiffraction->ProjFragment(pDef->GetPDGEncoding(),
92 thePrimary.Get4Momentum(),Z, A-Z);
93
94 if ( result && result->size() > 0 )
95 {
96 ktv=new G4KineticTrackVector();
97 std::vector<G4QHadron *>::iterator i;
98 #ifdef debug_scatter
99 G4int count(0);
100 #endif
101 for (i=result->begin(); i!=result->end(); i++)
102 {
103 G4int pdgCode=(*i)->GetPDGCode();
104 G4ParticleDefinition * secDef(0);
105 if ( pdgCode < 80000000) { // check if ion; should be 90Million, but 89Million is possible
107 }else {
108 G4int qN = pdgCode % 1000;
109 G4int qZ = (pdgCode/1000) %1000;
110 if ( qZ < 500 && qN < 500){ // protect for delta being coded as/in nucleus
111 secDef = G4ParticleTable::GetParticleTable()->GetIon(qZ,qN+qZ, 0);
112 if ( ! secDef )
113 { // exceptions to the rule!
114 if ( pdgCode == 90000001 ) secDef= G4Neutron::Neutron();
115 if ( pdgCode == 91000000 ) secDef= G4Lambda::Lambda();
116 }
117 }
118 }
119
120 if ( secDef ){
121
122 G4ThreeVector pos=(*i)->GetPosition(); //GetPosition returns const &
123 G4LorentzVector mom=(*i)->Get4Momentum();
124
125 G4KineticTrack * sPrim=new G4KineticTrack(secDef,
126 (*i)->GetFormationTime(), pos, mom);
127 ktv->push_back(sPrim);
128
129 #ifdef debug_scatter
130 G4cout << "G4ProjectileDiffractive sec # " << ++count << ", "
131 << "ChipsPDGCode=" << pdgCode << ", "
132 << secDef->GetParticleName()
133 << ", time(ns) " << (*i)->GetFormationTime()/ns
134 << ", pos " << pos
135 << ", 4mom " << (*i)->Get4Momentum()
136 << G4endl;
137 #endif
138 } else {
139 #ifdef debug_scatter
140 G4cout << "G4ProjectileDiffractiveChannel: Particle with PDG code "<< pdgCode <<" does not converted!!!"<<G4endl;
141 #endif
142 tryAgain=true;
143 }
144 }
145 std::for_each(result->begin(), result->end(), DeleteQHadron());
146 delete result;
147 result=0;
148 }
149
150 if ( result ) delete result;
151 if ( tryAgain ) { // QDiffraction returned a "difficult" particle in nucleus
152 std::for_each(ktv->begin(), ktv->end(), DeleteKineticTrack());
153 delete ktv;
154 ktv=0;
155 }
156 } while ( ! ktv);
157
158 return ktv;
159}
std::vector< G4QHadron * > G4QHadronVector
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
G4LorentzVector Get4Momentum() const
static G4Lambda * Lambda()
Definition: G4Lambda.cc:108
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
const G4String & GetParticleName() const
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()
G4ParticleDefinition * GetIon(G4int atomicNumber, G4int atomicMass, G4double excitationEnergy)
G4QHadronVector * ProjFragment(G4int pPDG, G4LorentzVector p4M, G4int tgZ, G4int tgN)
#define ns
Definition: xmlparse.cc:597

Referenced by G4TheoFSGenerator::ApplyYourself().


The documentation for this class was generated from the following files: