Changeset 1228 for trunk/source/processes/hadronic/models/qmd
- Timestamp:
- Jan 8, 2010, 11:56:51 AM (16 years ago)
- Location:
- trunk/source/processes/hadronic/models/qmd
- Files:
-
- 4 edited
-
GNUmakefile (modified) (3 diffs)
-
History (modified) (1 diff)
-
include/G4QMDReaction.hh (modified) (2 diffs)
-
src/G4QMDReaction.cc (modified) (8 diffs)
Legend:
- Unmodified
- Added
- Removed
-
trunk/source/processes/hadronic/models/qmd/GNUmakefile
r1055 r1228 1 # $Id: GNUmakefile,v 1. 2 2008/11/22 06:30:46tkoi Exp $1 # $Id: GNUmakefile,v 1.3 2009/11/20 23:33:38 tkoi Exp $ 2 2 # ----------------------------------------------------------- 3 3 # GNUmakefile for hadronic library. Gabriele Cosmo, 18/9/96. … … 30 30 -I$(G4BASE)/processes/hadronic/models/binary_cascade/include \ 31 31 -I$(G4BASE)/processes/hadronic/models/im_r_matrix/include \ 32 -I$(G4BASE)/processes/hadronic/models/chiral_inv_phase_space/interface/include \33 -I$(G4BASE)/processes/hadronic/models/chiral_inv_phase_space/body/include \34 32 -I$(G4BASE)/processes/hadronic/models/parton_string/qgsm/include \ 35 33 -I$(G4BASE)/processes/hadronic/models/parton_string/diffraction/include \ … … 43 41 -I$(G4BASE)/processes/hadronic/models/de_excitation/handler/include \ 44 42 -I$(G4BASE)/processes/hadronic/models/de_excitation/management/include \ 45 -I$(G4BASE)/processes/hadronic/models/pre_equilibrium/exciton_model/include \46 43 -I$(G4BASE)/particles/management/include \ 47 44 -I$(G4BASE)/particles/leptons/include \ -
trunk/source/processes/hadronic/models/qmd/History
r1196 r1228 13 13 * Please list in reverse chronological order (last date on top) 14 14 --------------------------------------------------------------- 15 16 20-November-2009 Tatsumi Koi (hadr-qmd-V09-02-03) 17 - Delete obviously unnecessary dependency in GNUMakefile 18 GNUmakefile 19 - Enable meson incidence (BETA) 20 include/G4QMDReaction.hh 21 src/G4QMDReaction.cc 15 22 16 23 19-November-2009 Tatsumi Koi (hadr-qmd-V09-02-02) -
trunk/source/processes/hadronic/models/qmd/include/G4QMDReaction.hh
r1055 r1228 50 50 #include "G4IonsShenCrossSection.hh" 51 51 #include "G4GeneralSpaceNNCrossSection.hh" 52 #include "G4PiNuclearCrossSection.hh" 52 53 53 54 #include "G4HadronicInteraction.hh" … … 115 116 G4GeneralSpaceNNCrossSection* genspaXS; 116 117 118 G4PiNuclearCrossSection* piNucXS; 119 117 120 G4bool gem; 118 121 G4bool frag; -
trunk/source/processes/hadronic/models/qmd/src/G4QMDReaction.cc
r1196 r1228 52 52 shenXS = new G4IonsShenCrossSection(); 53 53 //genspaXS = new G4GeneralSpaceNNCrossSection(); 54 piNucXS = new G4PiNuclearCrossSection(); 54 55 meanField = new G4QMDMeanField(); 55 56 collision = new G4QMDCollision(); … … 107 108 108 109 //090331 109 G4double xs_0 = shenXS->GetCrossSection ( proj_dp , targ_ele , aTemp ); 110 111 G4VCrossSectionDataSet* theXS = shenXS; 112 113 if ( proj_pd->GetParticleType() == "meson" ) theXS = piNucXS; 114 115 G4double xs_0 = theXS->GetCrossSection ( proj_dp , targ_ele , aTemp ); 116 110 117 //G4double xs_0 = genspaXS->GetCrossSection ( proj_dp , targ_ele , aTemp ); 111 118 G4double bmax_0 = std::sqrt( xs_0 / pi ); … … 158 165 159 166 G4double plab = projectile.GetTotalMomentum()/GeV; 160 G4double elab = ( projectile.GetKineticEnergy() + proj_pd->GetPDGMass() + targ_pd->GetPDGMass() )/GeV;167 G4double elab = ( projectile.GetKineticEnergy() + proj_pd->GetPDGMass() + targ_pd->GetPDGMass() )/GeV; 161 168 162 169 calcOffSetOfCollision( b , proj_pd , targ_pd , plab , elab , bmax , boostCMtoNN ); … … 177 184 //proj->ShowParticipants(); 178 185 179 } 180 181 meanField->SetSystem ( proj ); 182 proj->SetTotalPotential( meanField->GetTotalPotential() ); 183 proj->CalEnergyAndAngularMomentumInCM(); 186 187 meanField->SetSystem ( proj ); 188 proj->SetTotalPotential( meanField->GetTotalPotential() ); 189 proj->CalEnergyAndAngularMomentumInCM(); 190 191 } 184 192 185 193 // Target … … 220 228 G4LorentzVector targ4pCM = CLHEP::boostOf ( targ4p , boostToReac ); 221 229 222 223 230 // Projectile 224 231 if ( proj != NULL ) … … 251 258 // projectile is particle 252 259 253 G4int i = targ->GetTotalNumberOfParticipant(); 260 // avoid multiple set in "elastic" loop 261 if ( system->GetTotalNumberOfParticipant() == targ->GetTotalNumberOfParticipant() ) 262 { 263 264 G4int i = targ->GetTotalNumberOfParticipant(); 254 265 255 G4ThreeVector p0( 0 ); 256 G4ThreeVector r0( 0 ); 257 258 259 G4ThreeVector p ( p0.x() + coulomb_collision_px_proj 260 , p0.y() 261 , p0.z() * coulomb_collision_gamma_proj + coulomb_collision_pz_proj ); 262 263 G4ThreeVector r ( r0.x() + coulomb_collision_rx_proj 264 , r0.y() 265 , r0.z() / coulomb_collision_gamma_proj + coulomb_collision_rz_proj ); 266 267 system->SetParticipant( new G4QMDParticipant( (G4ParticleDefinition*)projectile.GetDefinition() , p , r ) ); 268 system->GetParticipant ( i )->SetProjectile(); 269 } 266 G4ThreeVector p0( 0 ); 267 G4ThreeVector r0( 0 ); 268 269 G4ThreeVector p ( p0.x() + coulomb_collision_px_proj 270 , p0.y() 271 , p0.z() * coulomb_collision_gamma_proj + coulomb_collision_pz_proj ); 272 273 G4ThreeVector r ( r0.x() + coulomb_collision_rx_proj 274 , r0.y() 275 , r0.z() / coulomb_collision_gamma_proj + coulomb_collision_rz_proj ); 276 277 system->SetParticipant( new G4QMDParticipant( (G4ParticleDefinition*)projectile.GetDefinition() , p , r ) ); 278 // This is not important becase only 1 projectile particle. 279 system->GetParticipant ( i )->SetProjectile(); 280 } 281 282 } 283 //system->ShowParticipants(); 270 284 271 285 delete targ; 272 286 delete proj; 273 274 287 275 288 meanField->SetSystem ( system ); … … 611 624 G4double ptot , G4double etot , G4double bmax , G4ThreeVector boostToCM ) 612 625 { 626 613 627 G4double mass_proj = pd_proj->GetPDGMass()/GeV; 614 628 G4double mass_targ = pd_targ->GetPDGMass()/GeV; … … 623 637 G4double eccm = stot - ( mass_proj + mass_targ ); 624 638 625 G4int zp = pd_proj->GetAtomicNumber(); 626 G4int ap = pd_proj->GetAtomicMass(); 639 G4int zp = 1; 640 G4int ap = 1; 641 if ( pd_proj->GetParticleType() == "nucleus" ) 642 { 643 zp = pd_proj->GetAtomicNumber(); 644 ap = pd_proj->GetAtomicMass(); 645 } 646 else 647 { 648 // proton, neutron, mesons 649 zp = int ( pd_proj->GetPDGCharge()/eplus + 0.5 ); 650 // ap = 1; 651 } 652 653 627 654 G4int zt = pd_targ->GetAtomicNumber(); 628 655 G4int at = pd_targ->GetAtomicMass(); 656 629 657 630 658 //G4double rmax0 = 8.0; // T.K dicide parameter value // for low energy
Note:
See TracChangeset
for help on using the changeset viewer.
