source: trunk/source/processes/hadronic/models/theo_high_energy/src/G4ProjectileDiffractiveChannel.cc@ 1199

Last change on this file since 1199 was 1196, checked in by garnier, 16 years ago

update CVS release candidate geant4.9.3.01

File size: 5.7 KB
Line 
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// $Id: G4ProjectileDiffractiveChannel.cc,v 1.2 2007/11/15 16:07:42 gunter Exp $
28// GEANT4 tag $Name: geant4-09-03-cand-01 $
29//
30
31// Author : Gunter Folger Nov 2007
32// Class Description
33// Final state production model for theoretical models of hadron inelastic
34// projectile diffractive scattering in geant4;
35// Class Description - End
36//
37// Modified:
38
39#include "G4ProjectileDiffractiveChannel.hh"
40
41#include "G4HadTmpUtil.hh" //lrint
42#include "G4QHadronVector.hh"
43#include "G4ParticleTable.hh"
44#include "G4IonTable.hh"
45#include "G4Lambda.hh"
46
47//#define debug_getFraction 1
48//#define debug_scatter 1
49
50G4ProjectileDiffractiveChannel::G4ProjectileDiffractiveChannel()
51{
52 theQDiffraction=G4QDiffractionRatio::GetPointer();
53}
54
55G4ProjectileDiffractiveChannel::~G4ProjectileDiffractiveChannel()
56{}
57
58G4double G4ProjectileDiffractiveChannel::GetFraction(G4Nucleus &theNucleus,
59 const G4DynamicParticle & thePrimary)
60{
61 G4double ratio;
62 ratio=theQDiffraction->GetRatio(thePrimary.GetTotalMomentum(),
63 thePrimary.GetDefinition()->GetPDGEncoding(),
64 G4lrint(theNucleus.GetZ()),
65 G4lrint(theNucleus.GetN()-theNucleus.GetZ()));
66 #ifdef debug_getFraction
67 G4cout << "G4ProjectileDiffractiveChannel::ratio " << ratio << G4endl;
68 #endif
69
70 return ratio;
71}
72
73G4KineticTrackVector * G4ProjectileDiffractiveChannel::Scatter(G4Nucleus &theNucleus,
74 const G4DynamicParticle & thePrimary)
75{
76
77
78 G4int A=G4lrint(theNucleus.GetN());
79 G4int Z=G4lrint(theNucleus.GetZ());
80
81 G4ParticleDefinition * pDef=thePrimary.GetDefinition();
82 #ifdef debug_scatter
83 G4cout << " ProjectileDiffractive: A, Z, proj-pdg" <<" "<< A <<" "<<Z
84 << " "<<" " << pDef->GetParticleName()<< G4endl;
85 #endif
86
87 G4QHadronVector * result(0);
88 G4KineticTrackVector * ktv(0);
89 G4bool tryAgain;
90
91 do {
92 tryAgain = false;
93 result=theQDiffraction->ProjFragment(pDef->GetPDGEncoding(),
94 thePrimary.Get4Momentum(),Z, A-Z);
95
96 if ( result && result->size() > 0 )
97 {
98 ktv=new G4KineticTrackVector();
99 std::vector<G4QHadron *>::iterator i;
100 #ifdef debug_scatter
101 G4int count(0);
102 #endif
103 for (i=result->begin(); i!=result->end(); i++)
104 {
105 G4int pdgCode=(*i)->GetPDGCode();
106 G4ParticleDefinition * secDef(0);
107 if ( pdgCode < 80000000) { // check if ion; should be 90Million, but 89Million is possible
108 secDef= G4ParticleTable::GetParticleTable()->FindParticle(pdgCode);
109 }else {
110 G4int N = pdgCode % 1000;
111 G4int Z = (pdgCode/1000) %1000;
112 if ( Z < 500 && N < 500){ // protect for delta being coded as/in nucleus
113 secDef = G4ParticleTable::GetParticleTable()->GetIon(Z,N+Z, 0);
114 if ( ! secDef )
115 { // exceptions to the rule!
116 if ( pdgCode == 90000001 ) secDef= G4Neutron::Neutron();
117 if ( pdgCode == 91000000 ) secDef= G4Lambda::Lambda();
118 }
119 }
120 }
121
122 if ( secDef ){
123
124 G4ThreeVector pos=(*i)->GetPosition(); //GetPosition returns const &
125 G4LorentzVector mom=(*i)->Get4Momentum();
126
127 G4KineticTrack * sPrim=new G4KineticTrack(secDef,
128 (*i)->GetFormationTime(), pos, mom);
129 ktv->push_back(sPrim);
130
131 #ifdef debug_scatter
132 G4cout << "G4ProjectileDiffractive sec # " << ++count << ", "
133 << "ChipsPDGCode=" << pdgCode << ", "
134 << secDef->GetParticleName()
135 << ", time(ns) " << (*i)->GetFormationTime()/ns
136 << ", pos " << pos
137 << ", 4mom " << (*i)->Get4Momentum()
138 << G4endl;
139 #endif
140 } else {
141 #ifdef debug_scatter
142 G4cout << "G4ProjectileDiffractiveChannel: Particle with PDG code "<< pdgCode <<" does not converted!!!"<<G4endl;
143 #endif
144 tryAgain=true;
145 }
146 }
147 std::for_each(result->begin(), result->end(), DeleteQHadron());
148 delete result;
149 result=0;
150 }
151
152 if ( result ) delete result;
153 if ( tryAgain ) { // QDiffraction returned a "difficult" particle in nucleus
154 std::for_each(ktv->begin(), ktv->end(), DeleteKineticTrack());
155 delete ktv;
156 ktv=0;
157 }
158 } while ( ! ktv);
159
160 return ktv;
161}
Note: See TracBrowser for help on using the repository browser.