source: trunk/source/processes/hadronic/models/incl/src/G4InclInput.cc@ 1357

Last change on this file since 1357 was 1350, checked in by garnier, 15 years ago

update to last version 4.9.4

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#include "G4InclInput.hh"
27
28G4InclInput::G4InclInput(const G4HadProjectile &aTrack, const G4Nucleus &theNucleus, G4bool inverseKinematics = false) {
29 usingInverseKinematics = inverseKinematics;
30
31 fTargetA = theNucleus.GetA_asInt(); // Target mass number
32 fTargetZ = theNucleus.GetZ_asInt(); // Target charge number
33 fBulletType = getBulletType(aTrack.GetDefinition()); // Projectile type (INCL particle code)
34 fBulletE = aTrack.GetKineticEnergy() / MeV; // Projectile energy (total, in MeV)
35 fTimeScale = 1.0; // Time scaling
36 fNuclearPotential = 45.0; // Nuclear potential
37 setExtendedProjectileInfo(aTrack.GetDefinition());
38 icoup = 0;
39 breakupThreshold = 10;
40 fMinNeutronEnergy = 0.0;
41 fMinProtonE = 0.0;
42}
43
44G4InclInput::~G4InclInput() {}
45
46void G4InclInput::printInfo() {
47 G4cout <<"Target: A = " << targetA() << " Z = " << targetZ() << G4endl;
48 G4cout <<"Projectile: type = " << bulletType() << " energy = " << bulletE() << G4endl;
49}
50
51void G4InclInput::printProjectileTargetInfo(const G4HadProjectile &aTrack, const G4Nucleus &theNucleus) {
52 G4cout <<"Projectile = " << aTrack.GetDefinition()->GetParticleName() << G4endl;
53 G4cout <<" four-momentum: " << aTrack.Get4Momentum() << G4endl;
54 G4cout <<"Energy = " << aTrack.GetKineticEnergy() / MeV << G4endl;
55 G4cout <<"Target A = " << theNucleus.GetA_asInt() << " Z = " << theNucleus.GetZ_asInt() << G4endl;
56}
57
58
59G4bool G4InclInput::canUseInverseKinematics(const G4HadProjectile &aTrack, const G4Nucleus &theNucleus) {
60 G4int targetA = theNucleus.GetA_asInt();
61 const G4ParticleDefinition *projectileDef = aTrack.GetDefinition();
62 G4int projectileA = projectileDef->GetAtomicMass();
63 // G4int projectileZ = projectileDef->GetAtomicNumber();
64 if(targetA > 0 && targetA < 18 && (projectileDef != G4Proton::Proton() &&
65 projectileDef != G4Neutron::Neutron() &&
66 projectileDef != G4PionPlus::PionPlus() &&
67 projectileDef != G4PionZero::PionZero() &&
68 projectileDef != G4PionMinus::PionMinus()) &&
69 projectileA > 1) {
70 return true;
71 } else {
72 return false;
73 }
74}
75
76void G4InclInput::setExtendedProjectileInfo(const G4ParticleDefinition *pd) {
77 if(getBulletType(pd) == -666) {
78 theExtendedProjectileA = pd->GetAtomicMass();
79 theExtendedProjectileZ = pd->GetAtomicNumber();
80 isExtended = true;
81 } else {
82 isExtended = false;
83 }
84}
85
86G4int G4InclInput::getBulletType(const G4ParticleDefinition *pd) {
87 // G4ParticleTable *pt = G4ParticleTable::GetParticleTable();
88
89 if(pd == G4Proton::Proton()) {
90 return 1;
91 } else if(pd == G4Neutron::Neutron()) {
92 return 2;
93 } else if(pd == G4PionPlus::PionPlus()) {
94 return 3;
95 } else if(pd == G4PionMinus::PionMinus()) {
96 return 5;
97 } else if(pd == G4PionZero::PionZero()) {
98 return 4;
99 } else if(pd == G4Deuteron::Deuteron()) {
100 return 6;
101 } else if(pd == G4Triton::Triton()) {
102 return 7;
103 } else if(pd == G4He3::He3()) {
104 return 8;
105 } else if(pd == G4Alpha::Alpha()) {
106 return 9;
107 // } else if(pd == pt->GetIon(6, 12, 0.0)) { // C12 special case. This should be phased-out in favor of "extended projectile"
108 // return -12;
109 } else { // Is this extended projectile?
110 G4int A = pd->GetAtomicMass();
111 G4int Z = pd->GetAtomicNumber();
112 if(A > 4 && A <= 16 && Z > 2 && Z <= 8) { // Ions from Lithium to Oxygen
113 return -666; // Code of an extended projectile
114 }
115 }
116 G4cout <<"Error! Projectile " << pd->GetParticleName() << " not defined!" << G4endl;
117 return 0;
118}
119
120G4ParticleDefinition* G4InclInput::getParticleDefinition(G4int inclParticleCode) {
121 switch(inclParticleCode) {
122 case 1:
123 return G4Proton::ProtonDefinition();
124 break;
125 case 2:
126 return G4Neutron::NeutronDefinition();
127 break;
128 case 3:
129 return G4PionPlus::PionPlusDefinition();
130 break;
131 case 4:
132 return G4PionMinus::PionMinusDefinition();
133 break;
134 case 5:
135 return G4PionZero::PionZeroDefinition();
136 break;
137 case 6:
138 return G4Deuteron::DeuteronDefinition();
139 break;
140 case 7:
141 return G4Triton::Triton();
142 break;
143 case 8:
144 return G4He3::He3Definition();
145 break;
146 case 9:
147 return G4Alpha::AlphaDefinition();
148 break;
149 }
150 return 0;
151}
Note: See TracBrowser for help on using the repository browser.