Geant4 Cross Reference |
1 // 1 // 2 // ******************************************* 2 // ******************************************************************** 3 // * License and Disclaimer 3 // * License and Disclaimer * 4 // * 4 // * * 5 // * The Geant4 software is copyright of th 5 // * The Geant4 software is copyright of the Copyright Holders of * 6 // * the Geant4 Collaboration. It is provided 6 // * the Geant4 Collaboration. It is provided under the terms and * 7 // * conditions of the Geant4 Software License 7 // * conditions of the Geant4 Software License, included in the file * 8 // * LICENSE and available at http://cern.ch/ 8 // * LICENSE and available at http://cern.ch/geant4/license . These * 9 // * include a list of copyright holders. 9 // * include a list of copyright holders. * 10 // * 10 // * * 11 // * Neither the authors of this software syst 11 // * Neither the authors of this software system, nor their employing * 12 // * institutes,nor the agencies providing fin 12 // * institutes,nor the agencies providing financial support for this * 13 // * work make any representation or warran 13 // * work make any representation or warranty, express or implied, * 14 // * regarding this software system or assum 14 // * regarding this software system or assume any liability for its * 15 // * use. Please see the license in the file 15 // * use. Please see the license in the file LICENSE and URL above * 16 // * for the full disclaimer and the limitatio 16 // * for the full disclaimer and the limitation of liability. * 17 // * 17 // * * 18 // * This code implementation is the result 18 // * This code implementation is the result of the scientific and * 19 // * technical work of the GEANT4 collaboratio 19 // * technical work of the GEANT4 collaboration. * 20 // * By using, copying, modifying or distri 20 // * By using, copying, modifying or distributing the software (or * 21 // * any work based on the software) you ag 21 // * any work based on the software) you agree to acknowledge its * 22 // * use in resulting scientific publicati 22 // * use in resulting scientific publications, and indicate your * 23 // * acceptance of all terms of the Geant4 Sof 23 // * acceptance of all terms of the Geant4 Software license. * 24 // ******************************************* 24 // ******************************************************************** 25 // 25 // 26 26 27 // G4 Low energy model: n-p scattering 27 // G4 Low energy model: n-p scattering 28 // F.W. Jones, L.G. Greeniaus, H.P. Wellisch 28 // F.W. Jones, L.G. Greeniaus, H.P. Wellisch 29 29 30 // 11-OCT-2007 F.W. Jones: removed erroneous c 30 // 11-OCT-2007 F.W. Jones: removed erroneous code for identity 31 // exchange of particles. 31 // exchange of particles. 32 // FWJ 27-AUG-2010: extended to 5 GeV by Tony 32 // FWJ 27-AUG-2010: extended to 5 GeV by Tony Kwan TRIUMF 33 // 06-Aug-15 A.Ribon migrating to G4Pow 33 // 06-Aug-15 A.Ribon migrating to G4Pow 34 34 35 #include "G4LEHadronProtonElastic.hh" 35 #include "G4LEHadronProtonElastic.hh" 36 #include "G4PhysicalConstants.hh" 36 #include "G4PhysicalConstants.hh" 37 #include "G4SystemOfUnits.hh" 37 #include "G4SystemOfUnits.hh" 38 #include "Randomize.hh" 38 #include "Randomize.hh" 39 #include "G4ios.hh" 39 #include "G4ios.hh" 40 #include "G4Pow.hh" 40 #include "G4Pow.hh" 41 #include "G4PhysicsModelCatalog.hh" 41 #include "G4PhysicsModelCatalog.hh" 42 42 43 43 44 G4LEHadronProtonElastic::G4LEHadronProtonElast 44 G4LEHadronProtonElastic::G4LEHadronProtonElastic(): 45 G4HadronElastic("G4LEHadronProtonElastic") 45 G4HadronElastic("G4LEHadronProtonElastic") 46 { 46 { 47 secID = G4PhysicsModelCatalog::GetModelID( " 47 secID = G4PhysicsModelCatalog::GetModelID( "model_" + GetModelName() ); 48 SetMinEnergy(0.); 48 SetMinEnergy(0.); 49 SetMaxEnergy(20.*MeV); 49 SetMaxEnergy(20.*MeV); 50 } 50 } 51 51 52 G4LEHadronProtonElastic::~G4LEHadronProtonElas 52 G4LEHadronProtonElastic::~G4LEHadronProtonElastic() 53 { 53 { 54 theParticleChange.Clear(); 54 theParticleChange.Clear(); 55 } 55 } 56 56 57 G4HadFinalState* 57 G4HadFinalState* 58 G4LEHadronProtonElastic::ApplyYourself(const G 58 G4LEHadronProtonElastic::ApplyYourself(const G4HadProjectile& aTrack, 59 G4Nucleus& targetN 59 G4Nucleus& targetNucleus) 60 { 60 { 61 theParticleChange.Clear(); 61 theParticleChange.Clear(); 62 const G4HadProjectile* aParticle = &aTrack 62 const G4HadProjectile* aParticle = &aTrack; 63 63 64 G4double P = aParticle->GetTotalMomentum() 64 G4double P = aParticle->GetTotalMomentum(); 65 G4double Px = aParticle->Get4Momentum().x( 65 G4double Px = aParticle->Get4Momentum().x(); 66 G4double Py = aParticle->Get4Momentum().y( 66 G4double Py = aParticle->Get4Momentum().y(); 67 G4double Pz = aParticle->Get4Momentum().z( 67 G4double Pz = aParticle->Get4Momentum().z(); 68 G4double ek = aParticle->GetKineticEnergy( 68 G4double ek = aParticle->GetKineticEnergy(); 69 G4ThreeVector theInitial = aParticle->Get4 69 G4ThreeVector theInitial = aParticle->Get4Momentum().vect(); 70 70 71 if (verboseLevel > 1) 71 if (verboseLevel > 1) 72 { 72 { 73 G4double E = aParticle->GetTotalEnergy() 73 G4double E = aParticle->GetTotalEnergy(); 74 G4double E0 = aParticle->GetDefinition() 74 G4double E0 = aParticle->GetDefinition()->GetPDGMass(); 75 G4double Q = aParticle->GetDefinition()- 75 G4double Q = aParticle->GetDefinition()->GetPDGCharge(); 76 G4int A = targetNucleus.GetA_asInt(); 76 G4int A = targetNucleus.GetA_asInt(); 77 G4int Z = targetNucleus.GetZ_asInt(); 77 G4int Z = targetNucleus.GetZ_asInt(); 78 G4cout << "G4LEHadronProtonElastic:Apply 78 G4cout << "G4LEHadronProtonElastic:ApplyYourself: incident particle: " 79 << aParticle->GetDefinition()->Ge 79 << aParticle->GetDefinition()->GetParticleName() << G4endl; 80 G4cout << "P = " << P/GeV << " GeV/c" 80 G4cout << "P = " << P/GeV << " GeV/c" 81 << ", Px = " << Px/GeV << " GeV/c 81 << ", Px = " << Px/GeV << " GeV/c" 82 << ", Py = " << Py/GeV << " GeV/c 82 << ", Py = " << Py/GeV << " GeV/c" 83 << ", Pz = " << Pz/GeV << " GeV/c 83 << ", Pz = " << Pz/GeV << " GeV/c" << G4endl; 84 G4cout << "E = " << E/GeV << " GeV" 84 G4cout << "E = " << E/GeV << " GeV" 85 << ", kinetic energy = " << ek/Ge 85 << ", kinetic energy = " << ek/GeV << " GeV" 86 << ", mass = " << E0/GeV << " GeV 86 << ", mass = " << E0/GeV << " GeV" 87 << ", charge = " << Q << G4endl; 87 << ", charge = " << Q << G4endl; 88 G4cout << "G4LEHadronProtonElastic:Apply 88 G4cout << "G4LEHadronProtonElastic:ApplyYourself: material:" << G4endl; 89 G4cout << "A = " << A 89 G4cout << "A = " << A 90 << ", Z = " << Z 90 << ", Z = " << Z 91 << ", atomic mass " 91 << ", atomic mass " 92 << G4Proton::Proton()->GetPDGMas 92 << G4Proton::Proton()->GetPDGMass()/GeV << "GeV" 93 << G4endl; 93 << G4endl; 94 // 94 // 95 // GHEISHA ADD operation to get total en 95 // GHEISHA ADD operation to get total energy, mass, charge 96 // 96 // 97 E += proton_mass_c2; 97 E += proton_mass_c2; 98 G4double E02 = E*E - P*P; 98 G4double E02 = E*E - P*P; 99 E0 = std::sqrt(std::abs(E02)); 99 E0 = std::sqrt(std::abs(E02)); 100 if (E02 < 0)E0 *= -1; 100 if (E02 < 0)E0 *= -1; 101 Q += Z; 101 Q += Z; 102 G4cout << "G4LEHadronProtonElastic:Apply 102 G4cout << "G4LEHadronProtonElastic:ApplyYourself: total:" << G4endl; 103 G4cout << "E = " << E/GeV << " GeV" 103 G4cout << "E = " << E/GeV << " GeV" 104 << ", mass = " << E0/GeV << " GeV 104 << ", mass = " << E0/GeV << " GeV" 105 << ", charge = " << Q << G4endl; 105 << ", charge = " << Q << G4endl; 106 } 106 } 107 107 108 G4double theta = (0.5)*pi/180.; 108 G4double theta = (0.5)*pi/180.; 109 109 110 // Get the target particle 110 // Get the target particle 111 111 112 G4DynamicParticle* targetParticle = target 112 G4DynamicParticle* targetParticle = targetNucleus.ReturnTargetParticle(); 113 113 114 G4double E1 = aParticle->GetTotalEnergy(); 114 G4double E1 = aParticle->GetTotalEnergy(); 115 G4double M1 = aParticle->GetDefinition()-> 115 G4double M1 = aParticle->GetDefinition()->GetPDGMass(); 116 G4double E2 = targetParticle->GetTotalEner 116 G4double E2 = targetParticle->GetTotalEnergy(); 117 G4double M2 = targetParticle->GetDefinitio 117 G4double M2 = targetParticle->GetDefinition()->GetPDGMass(); 118 G4double totalEnergy = E1 + E2; 118 G4double totalEnergy = E1 + E2; 119 G4double pseudoMass = std::sqrt(totalEnerg 119 G4double pseudoMass = std::sqrt(totalEnergy*totalEnergy - P*P); 120 120 121 // Transform into centre of mass system 121 // Transform into centre of mass system 122 122 123 G4double px = (M2/pseudoMass)*Px; 123 G4double px = (M2/pseudoMass)*Px; 124 G4double py = (M2/pseudoMass)*Py; 124 G4double py = (M2/pseudoMass)*Py; 125 G4double pz = (M2/pseudoMass)*Pz; 125 G4double pz = (M2/pseudoMass)*Pz; 126 G4double p = std::sqrt(px*px + py*py + pz* 126 G4double p = std::sqrt(px*px + py*py + pz*pz); 127 127 128 if (verboseLevel > 1) { 128 if (verboseLevel > 1) { 129 G4cout << " E1, M1 (GeV) " << E1/GeV << 129 G4cout << " E1, M1 (GeV) " << E1/GeV << " " << M1/GeV << G4endl; 130 G4cout << " E2, M2 (GeV) " << E2/GeV << 130 G4cout << " E2, M2 (GeV) " << E2/GeV << " " << M2/GeV << G4endl; 131 G4cout << " particle 1 momentum in CM 131 G4cout << " particle 1 momentum in CM " << px/GeV << " " << py/GeV << " " 132 << pz/GeV << " " << p/GeV << G4endl 132 << pz/GeV << " " << p/GeV << G4endl; 133 } 133 } 134 134 135 // First scatter w.r.t. Z axis 135 // First scatter w.r.t. Z axis 136 G4double phi = G4UniformRand()*twopi; 136 G4double phi = G4UniformRand()*twopi; 137 G4double pxnew = p*std::sin(theta)*std::co 137 G4double pxnew = p*std::sin(theta)*std::cos(phi); 138 G4double pynew = p*std::sin(theta)*std::si 138 G4double pynew = p*std::sin(theta)*std::sin(phi); 139 G4double pznew = p*std::cos(theta); 139 G4double pznew = p*std::cos(theta); 140 140 141 // Rotate according to the direction of th 141 // Rotate according to the direction of the incident particle 142 if (px*px + py*py > 0) 142 if (px*px + py*py > 0) 143 { 143 { 144 G4double cost, sint, ph, cosp, sinp; 144 G4double cost, sint, ph, cosp, sinp; 145 cost = pz/p; 145 cost = pz/p; 146 sint = (std::sqrt(std::fabs((1-cost)*(1+ 146 sint = (std::sqrt(std::fabs((1-cost)*(1+cost))) 147 + std::sqrt(px*px+py*py)/p)/2; 147 + std::sqrt(px*px+py*py)/p)/2; 148 py < 0 ? ph = 3*halfpi : ph = halfpi; 148 py < 0 ? ph = 3*halfpi : ph = halfpi; 149 if (std::abs(px) > 0.000001*GeV) ph = st 149 if (std::abs(px) > 0.000001*GeV) ph = std::atan2(py,px); 150 cosp = std::cos(ph); 150 cosp = std::cos(ph); 151 sinp = std::sin(ph); 151 sinp = std::sin(ph); 152 px = (cost*cosp*pxnew - sinp*pynew + sin 152 px = (cost*cosp*pxnew - sinp*pynew + sint*cosp*pznew); 153 py = (cost*sinp*pxnew + cosp*pynew + sin 153 py = (cost*sinp*pxnew + cosp*pynew + sint*sinp*pznew); 154 pz = (-sint*pxnew + cos 154 pz = (-sint*pxnew + cost*pznew); 155 } 155 } 156 else { 156 else { 157 px = pxnew; 157 px = pxnew; 158 py = pynew; 158 py = pynew; 159 pz = pznew; 159 pz = pznew; 160 } 160 } 161 161 162 if (verboseLevel > 1) { 162 if (verboseLevel > 1) { 163 G4cout << " AFTER SCATTER..." << G4endl 163 G4cout << " AFTER SCATTER..." << G4endl; 164 G4cout << " particle 1 momentum in CM " 164 G4cout << " particle 1 momentum in CM " << px/GeV 165 << " " << py/GeV << " " << pz/GeV 165 << " " << py/GeV << " " << pz/GeV << " " << p/GeV 166 << G4endl; 166 << G4endl; 167 } 167 } 168 168 169 // Transform to lab system 169 // Transform to lab system 170 170 171 G4double E1pM2 = E1 + M2; 171 G4double E1pM2 = E1 + M2; 172 G4double betaCM = P/E1pM2; 172 G4double betaCM = P/E1pM2; 173 G4double betaCMx = Px/E1pM2; 173 G4double betaCMx = Px/E1pM2; 174 G4double betaCMy = Py/E1pM2; 174 G4double betaCMy = Py/E1pM2; 175 G4double betaCMz = Pz/E1pM2; 175 G4double betaCMz = Pz/E1pM2; 176 G4double gammaCM = E1pM2/std::sqrt(E1pM2*E 176 G4double gammaCM = E1pM2/std::sqrt(E1pM2*E1pM2 - P*P); 177 177 178 if (verboseLevel > 1) { 178 if (verboseLevel > 1) { 179 G4cout << " betaCM " << betaCMx << " " 179 G4cout << " betaCM " << betaCMx << " " << betaCMy << " " 180 << betaCMz << " " << betaCM << G4 180 << betaCMz << " " << betaCM << G4endl; 181 G4cout << " gammaCM " << gammaCM << G4e 181 G4cout << " gammaCM " << gammaCM << G4endl; 182 } 182 } 183 183 184 // Now following GLOREN... 184 // Now following GLOREN... 185 185 186 G4double BETA[5], PA[5], PB[5]; 186 G4double BETA[5], PA[5], PB[5]; 187 BETA[1] = -betaCMx; 187 BETA[1] = -betaCMx; 188 BETA[2] = -betaCMy; 188 BETA[2] = -betaCMy; 189 BETA[3] = -betaCMz; 189 BETA[3] = -betaCMz; 190 BETA[4] = gammaCM; 190 BETA[4] = gammaCM; 191 191 192 //The incident particle... 192 //The incident particle... 193 193 194 PA[1] = px; 194 PA[1] = px; 195 PA[2] = py; 195 PA[2] = py; 196 PA[3] = pz; 196 PA[3] = pz; 197 PA[4] = std::sqrt(M1*M1 + p*p); 197 PA[4] = std::sqrt(M1*M1 + p*p); 198 198 199 G4double BETPA = BETA[1]*PA[1] + BETA[2]* 199 G4double BETPA = BETA[1]*PA[1] + BETA[2]*PA[2] + BETA[3]*PA[3]; 200 G4double BPGAM = (BETPA * BETA[4]/(BETA[4 200 G4double BPGAM = (BETPA * BETA[4]/(BETA[4] + 1.) - PA[4]) * BETA[4]; 201 201 202 PB[1] = PA[1] + BPGAM * BETA[1]; 202 PB[1] = PA[1] + BPGAM * BETA[1]; 203 PB[2] = PA[2] + BPGAM * BETA[2]; 203 PB[2] = PA[2] + BPGAM * BETA[2]; 204 PB[3] = PA[3] + BPGAM * BETA[3]; 204 PB[3] = PA[3] + BPGAM * BETA[3]; 205 PB[4] = (PA[4] - BETPA) * BETA[4]; 205 PB[4] = (PA[4] - BETPA) * BETA[4]; 206 206 207 G4DynamicParticle* newP = new G4DynamicPar 207 G4DynamicParticle* newP = new G4DynamicParticle; 208 newP->SetDefinition(aParticle->GetDefiniti 208 newP->SetDefinition(aParticle->GetDefinition()); 209 newP->SetMomentum(G4ThreeVector(PB[1], PB[ 209 newP->SetMomentum(G4ThreeVector(PB[1], PB[2], PB[3])); 210 210 211 //The target particle... 211 //The target particle... 212 212 213 PA[1] = -px; 213 PA[1] = -px; 214 PA[2] = -py; 214 PA[2] = -py; 215 PA[3] = -pz; 215 PA[3] = -pz; 216 PA[4] = std::sqrt(M2*M2 + p*p); 216 PA[4] = std::sqrt(M2*M2 + p*p); 217 217 218 BETPA = BETA[1]*PA[1] + BETA[2]*PA[2] + B 218 BETPA = BETA[1]*PA[1] + BETA[2]*PA[2] + BETA[3]*PA[3]; 219 BPGAM = (BETPA * BETA[4]/(BETA[4] + 1.) - 219 BPGAM = (BETPA * BETA[4]/(BETA[4] + 1.) - PA[4]) * BETA[4]; 220 220 221 PB[1] = PA[1] + BPGAM * BETA[1]; 221 PB[1] = PA[1] + BPGAM * BETA[1]; 222 PB[2] = PA[2] + BPGAM * BETA[2]; 222 PB[2] = PA[2] + BPGAM * BETA[2]; 223 PB[3] = PA[3] + BPGAM * BETA[3]; 223 PB[3] = PA[3] + BPGAM * BETA[3]; 224 PB[4] = (PA[4] - BETPA) * BETA[4]; 224 PB[4] = (PA[4] - BETPA) * BETA[4]; 225 225 226 targetParticle->SetMomentum(G4ThreeVector( 226 targetParticle->SetMomentum(G4ThreeVector(PB[1], PB[2], PB[3])); 227 227 228 if (verboseLevel > 1) { 228 if (verboseLevel > 1) { 229 G4cout << " particle 1 momentum in LAB 229 G4cout << " particle 1 momentum in LAB " 230 << newP->GetMomentum()*(1./GeV) 230 << newP->GetMomentum()*(1./GeV) 231 << " " << newP->GetTotalMomentum()/ 231 << " " << newP->GetTotalMomentum()/GeV << G4endl; 232 G4cout << " particle 2 momentum in LAB 232 G4cout << " particle 2 momentum in LAB " 233 << targetParticle->GetMomentum()*(1 233 << targetParticle->GetMomentum()*(1./GeV) 234 << " " << targetParticle->GetTotalM 234 << " " << targetParticle->GetTotalMomentum()/GeV << G4endl; 235 G4cout << " TOTAL momentum in LAB " 235 G4cout << " TOTAL momentum in LAB " 236 << (newP->GetMomentum()+targetParti 236 << (newP->GetMomentum()+targetParticle->GetMomentum())*(1./GeV) 237 << " " 237 << " " 238 << (newP->GetMomentum()+targetParti 238 << (newP->GetMomentum()+targetParticle->GetMomentum()).mag()/GeV 239 << G4endl; 239 << G4endl; 240 } 240 } 241 241 242 theParticleChange.SetMomentumChange(newP-> 242 theParticleChange.SetMomentumChange(newP->GetMomentumDirection()); 243 theParticleChange.SetEnergyChange(newP->Ge 243 theParticleChange.SetEnergyChange(newP->GetKineticEnergy()); 244 delete newP; 244 delete newP; 245 theParticleChange.AddSecondary(targetParti 245 theParticleChange.AddSecondary(targetParticle, secID); 246 246 247 return &theParticleChange; 247 return &theParticleChange; 248 } 248 } 249 249 250 ////////////////////////////////////////////// 250 //////////////////////////////////////////////////////////////////// 251 // 251 // 252 // sample momentum transfer using Lab. momentu 252 // sample momentum transfer using Lab. momentum 253 253 254 G4double 254 G4double 255 G4LEHadronProtonElastic::SampleInvariantT(cons 255 G4LEHadronProtonElastic::SampleInvariantT(const G4ParticleDefinition* p, 256 G4double plab, G4int , G4int ) 256 G4double plab, G4int , G4int ) 257 { 257 { 258 G4double hMass = p->GetPDGMass(); 258 G4double hMass = p->GetPDGMass(); 259 G4double pCMS = 0.5*plab; 259 G4double pCMS = 0.5*plab; 260 // pCMS *= 50; 260 // pCMS *= 50; 261 G4double hEcms = std::sqrt(pCMS*pCMS+hMass*h 261 G4double hEcms = std::sqrt(pCMS*pCMS+hMass*hMass); 262 // G4double gamma = hEcms/hMass; 262 // G4double gamma = hEcms/hMass; 263 // gamma *= 15; 263 // gamma *= 15; 264 G4double beta = pCMS/hEcms; // std::sqrt(1- 264 G4double beta = pCMS/hEcms; // std::sqrt(1-1./gamma/gamma); // 265 // beta /= 0.8; // 0.95; // 1.0; // 1.1 // 0 265 // beta /= 0.8; // 0.95; // 1.0; // 1.1 // 0.5*pi; // pi; twopi; 266 G4double cosDipole = RandCosThetaDipPen(); 266 G4double cosDipole = RandCosThetaDipPen(); 267 G4double cosTheta = cosDipole + beta; 267 G4double cosTheta = cosDipole + beta; 268 cosTheta /= 1. + cosDipole*beta; 268 cosTheta /= 1. + cosDipole*beta; 269 G4double t = 2.*pCMS*pCMS*(1.-cosTheta); 269 G4double t = 2.*pCMS*pCMS*(1.-cosTheta); 270 return t; 270 return t; 271 271 272 } 272 } 273 273 274 ////////////////////////////////////////////// 274 /////////////////////////////////////////////////////////////// 275 // 275 // 276 // 1 + cos^2(theta) random distribution in the 276 // 1 + cos^2(theta) random distribution in the projectile rest frame, Penelope algorithm 277 277 278 G4double G4LEHadronProtonElastic::RandCosTheta 278 G4double G4LEHadronProtonElastic::RandCosThetaDipPen() 279 { 279 { 280 G4double x, cosTheta, signX, modX, power = 1 280 G4double x, cosTheta, signX, modX, power = 1./3.; 281 281 282 if( G4UniformRand() > 0.25) 282 if( G4UniformRand() > 0.25) 283 { 283 { 284 cosTheta = 2.*G4UniformRand()-1.; 284 cosTheta = 2.*G4UniformRand()-1.; 285 } 285 } 286 else 286 else 287 { 287 { 288 x = 2.*G4UniformRand()-1.; 288 x = 2.*G4UniformRand()-1.; 289 289 290 if ( x < 0. ) 290 if ( x < 0. ) 291 { 291 { 292 modX = -x; 292 modX = -x; 293 signX = -1.; 293 signX = -1.; 294 } 294 } 295 else 295 else 296 { 296 { 297 modX = x; 297 modX = x; 298 signX = 1.; 298 signX = 1.; 299 } 299 } 300 cosTheta = signX*G4Pow::GetInstance()->pow 300 cosTheta = signX*G4Pow::GetInstance()->powA(modX,power); 301 } 301 } 302 return cosTheta; 302 return cosTheta; 303 } 303 } 304 304 305 // end of file 305 // end of file 306 306