00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033 #define INCLXX_IN_GEANT4_MODE 1
00034
00035 #include "globals.hh"
00036
00037 #include <cmath>
00038
00039 #include "G4PhysicalConstants.hh"
00040 #include "G4SystemOfUnits.hh"
00041 #include "G4INCLXXInterface.hh"
00042 #include "G4GenericIon.hh"
00043 #include "G4INCLCascade.hh"
00044 #include "G4ReactionProductVector.hh"
00045 #include "G4ReactionProduct.hh"
00046 #include "G4INCLXXInterfaceStore.hh"
00047 #include "G4String.hh"
00048
00049 G4INCLXXInterface::G4INCLXXInterface(const G4String& nam) :
00050 G4VIntraNuclearTransportModel(nam),
00051 theINCLModel(NULL),
00052 theInterfaceStore(G4INCLXXInterfaceStore::GetInstance()),
00053 complainedAboutBackupModel(false)
00054 {
00055
00056 if(getenv("G4INCLXX_NO_DE_EXCITATION")) {
00057 G4String message = "de-excitation is completely disabled!";
00058 theInterfaceStore->EmitWarning(message);
00059 theExcitationHandler = 0;
00060 } else {
00061 theExcitationHandler = new G4ExcitationHandler;
00062 }
00063
00064 theBackupModel = new G4BinaryLightIonReaction;
00065 }
00066
00067 G4INCLXXInterface::~G4INCLXXInterface()
00068 {
00069 delete theBackupModel;
00070 delete theExcitationHandler;
00071 }
00072
00073 G4bool G4INCLXXInterface::AccurateProjectile(const G4HadProjectile &aTrack, const G4Nucleus &theNucleus) const {
00074
00075 const G4ParticleDefinition *projectileDef = aTrack.GetDefinition();
00076 if(projectileDef == G4Proton::Proton()
00077 || projectileDef == G4Neutron::Neutron()
00078 || projectileDef == G4PionPlus::PionPlus()
00079 || projectileDef == G4PionZero::PionZero()
00080 || projectileDef == G4PionMinus::PionMinus())
00081 return false;
00082
00083
00084 const G4int pA = projectileDef->GetAtomicMass();
00085 if(pA<=0) {
00086 std::stringstream ss;
00087 ss << "the model does not know how to handle a collision between a "
00088 << projectileDef->GetParticleName() << " projectile and a Z="
00089 << theNucleus.GetZ_asInt() << ", A=" << theNucleus.GetA_asInt();
00090 theInterfaceStore->EmitWarning(ss.str());
00091 return true;
00092 }
00093
00094
00095 const G4int tA = theNucleus.GetA_asInt();
00096 if(tA<=4 || pA<=4) {
00097 if(pA<tA)
00098 return false;
00099 else
00100 return true;
00101 }
00102
00103
00104
00105
00106
00107
00108 const G4int theMaxProjMassINCL = theInterfaceStore->GetMaxProjMassINCL();
00109 if(pA > theMaxProjMassINCL)
00110 return true;
00111 else if(tA > theMaxProjMassINCL)
00112 return false;
00113 else
00114
00115 return theInterfaceStore->GetAccurateProjectile();
00116 }
00117
00118 G4HadFinalState* G4INCLXXInterface::ApplyYourself(const G4HadProjectile& aTrack, G4Nucleus& theNucleus)
00119 {
00120
00121
00122 const G4int theMaxProjMassINCL = theInterfaceStore->GetMaxProjMassINCL();
00123 if(aTrack.GetDefinition()->GetAtomicMass() > theMaxProjMassINCL
00124 && theNucleus.GetA_asInt() > theMaxProjMassINCL) {
00125 if(!complainedAboutBackupModel) {
00126 complainedAboutBackupModel = true;
00127 std::stringstream ss;
00128 ss << "INCL++ refuses to handle reactions between nuclei with A>"
00129 << theMaxProjMassINCL
00130 << ". A backup model ("
00131 << theBackupModel->GetModelName()
00132 << ") will be used instead.";
00133 G4cout << "[INCL++] Warning: " << ss.str() << G4endl;
00134 }
00135 return theBackupModel->ApplyYourself(aTrack, theNucleus);
00136 }
00137
00138 const G4int maxTries = 200;
00139
00140
00141 const G4bool inverseKinematics = AccurateProjectile(aTrack, theNucleus);
00142
00143
00144 G4LorentzRotation *toInverseKinematics = NULL;
00145 G4LorentzRotation *toDirectKinematics = NULL;
00146 G4HadProjectile const *aProjectileTrack = &aTrack;
00147 G4Nucleus *theTargetNucleus = &theNucleus;
00148 if(inverseKinematics) {
00149 G4ParticleTable * const theParticleTable = G4ParticleTable::GetParticleTable();
00150 const G4int oldTargetA = theNucleus.GetA_asInt();
00151 const G4int oldTargetZ = theNucleus.GetZ_asInt();
00152 G4ParticleDefinition *oldTargetDef = theParticleTable->GetIon(oldTargetZ, oldTargetA, 0.0);
00153 const G4ParticleDefinition *oldProjectileDef = aTrack.GetDefinition();
00154
00155 if(oldProjectileDef != 0 && oldTargetDef != 0) {
00156 const G4int newTargetA = oldProjectileDef->GetAtomicMass();
00157 const G4int newTargetZ = oldProjectileDef->GetAtomicNumber();
00158
00159 if(newTargetA > 0 && newTargetZ > 0) {
00160
00161 theTargetNucleus = new G4Nucleus(newTargetA, newTargetZ);
00162 const G4double theProjectileMass = theParticleTable->GetIonTable()->GetIonMass(oldTargetZ, oldTargetA);
00163 toInverseKinematics = new G4LorentzRotation(aTrack.Get4Momentum().boostVector());
00164 G4LorentzVector theProjectile4Momentum(0.0, 0.0, 0.0, theProjectileMass);
00165 G4DynamicParticle swappedProjectileParticle(oldTargetDef, (*toInverseKinematics) * theProjectile4Momentum);
00166 aProjectileTrack = new G4HadProjectile(swappedProjectileParticle);
00167 } else {
00168 G4String message = "badly defined target after swapping. Falling back to normal (non-swapped) mode.";
00169 theInterfaceStore->EmitWarning(message);
00170 toInverseKinematics = new G4LorentzRotation;
00171 }
00172 } else {
00173 G4String message = "oldProjectileDef or oldTargetDef was null";
00174 theInterfaceStore->EmitWarning(message);
00175 toInverseKinematics = new G4LorentzRotation;
00176 }
00177 }
00178
00179
00180
00181
00182
00183 G4LorentzVector projectileMomentum = aProjectileTrack->Get4Momentum();
00184
00185
00186
00187
00188
00189
00190
00191
00192
00193 G4RotationMatrix toZ;
00194 toZ.rotateZ(-projectileMomentum.phi());
00195 toZ.rotateY(-projectileMomentum.theta());
00196 G4RotationMatrix toLabFrame3 = toZ.inverse();
00197 G4LorentzRotation toLabFrame(toLabFrame3);
00198
00199
00200
00201
00202
00203
00204
00205
00206 theResult.Clear();
00207 theResult.SetStatusChange(stopAndKill);
00208
00209 std::list<G4Fragment> remnants;
00210
00211 G4int nTries = 0;
00212
00213
00214
00215 G4bool eventIsOK = false;
00216 do {
00217 const G4INCL::ParticleSpecies theSpecies = toINCLParticleSpecies(*aProjectileTrack);
00218 const G4double kineticEnergy = toINCLKineticEnergy(*aProjectileTrack);
00219
00220
00221 theINCLModel = G4INCLXXInterfaceStore::GetInstance()->GetINCLModel();
00222
00223 if(theInterfaceStore->GetDumpInput()) {
00224 G4cout << theINCLModel->configToString() << G4endl;
00225 }
00226 const G4INCL::EventInfo eventInfo = theINCLModel->processEvent(theSpecies, kineticEnergy, theTargetNucleus->GetA_asInt(), theTargetNucleus->GetZ_asInt());
00227
00228 eventIsOK = !eventInfo.transparent;
00229 if(eventIsOK) {
00230
00231
00232
00233 if(inverseKinematics) {
00234 toDirectKinematics = new G4LorentzRotation(toInverseKinematics->inverse());
00235 }
00236
00237 for(G4int i = 0; i < eventInfo.nParticles; i++) {
00238 G4int A = eventInfo.A[i];
00239 G4int Z = eventInfo.Z[i];
00240
00241 G4double kinE = eventInfo.EKin[i];
00242 G4double px = eventInfo.px[i];
00243 G4double py = eventInfo.py[i];
00244 G4double pz = eventInfo.pz[i];
00245 G4DynamicParticle *p = toG4Particle(A, Z , kinE, px, py, pz);
00246 if(p != 0) {
00247 G4LorentzVector momentum = p->Get4Momentum();
00248
00249
00250 momentum *= toLabFrame;
00251
00252 if(inverseKinematics) {
00253 momentum *= *toDirectKinematics;
00254 momentum.setVect(-momentum.vect());
00255 }
00256
00257
00258 p->Set4Momentum(momentum);
00259 theResult.AddSecondary(p);
00260
00261 } else {
00262 G4String message = "the model produced a particle that couldn't be converted to Geant4 particle.";
00263 theInterfaceStore->EmitWarning(message);
00264 }
00265 }
00266
00267 for(G4int i = 0; i < eventInfo.nRemnants; i++) {
00268 const G4int A = eventInfo.ARem[i];
00269 const G4int Z = eventInfo.ZRem[i];
00270
00271 const G4double kinE = eventInfo.EKinRem[i];
00272 const G4double px = eventInfo.pxRem[i];
00273 const G4double py = eventInfo.pyRem[i];
00274 const G4double pz = eventInfo.pzRem[i];
00275 G4ThreeVector spin(
00276 eventInfo.jxRem[i]*hbar_Planck,
00277 eventInfo.jyRem[i]*hbar_Planck,
00278 eventInfo.jzRem[i]*hbar_Planck
00279 );
00280 const G4double excitationE = eventInfo.EStarRem[i];
00281 const G4double nuclearMass = G4NucleiProperties::GetNuclearMass(A, Z) + excitationE;
00282 const G4double scaling = remnant4MomentumScaling(nuclearMass,
00283 kinE,
00284 px, py, pz);
00285 G4LorentzVector fourMomentum(scaling * px, scaling * py, scaling * pz,
00286 nuclearMass + kinE);
00287 if(std::abs(scaling - 1.0) > 0.01) {
00288 std::stringstream ss;
00289 ss << "momentum scaling = " << scaling
00290 << "\n Lorentz vector = " << fourMomentum
00291 << ")\n A = " << A << ", Z = " << Z
00292 << "\n E* = " << excitationE << ", nuclearMass = " << nuclearMass
00293 << "\n remnant i=" << i << ", nRemnants=" << eventInfo.nRemnants
00294 << "\n Reaction was: " << aTrack.GetKineticEnergy()/MeV
00295 << "-MeV " << aTrack.GetDefinition()->GetParticleName() << " + "
00296 << G4ParticleTable::GetParticleTable()->GetIon(theNucleus.GetZ_asInt(), theNucleus.GetA_asInt(), 0.0)->GetParticleName()
00297 << ", in " << (inverseKinematics ? "inverse" : "direct") << " kinematics.";
00298 theInterfaceStore->EmitWarning(ss.str());
00299 }
00300
00301
00302 fourMomentum *= toLabFrame;
00303 spin *= toLabFrame3;
00304
00305 if(inverseKinematics) {
00306 fourMomentum *= *toDirectKinematics;
00307 fourMomentum.setVect(-fourMomentum.vect());
00308 }
00309
00310 G4Fragment remnant(A, Z, fourMomentum);
00311 remnant.SetAngularMomentum(spin);
00312 remnants.push_back(remnant);
00313 }
00314 }
00315 nTries++;
00316 } while(!eventIsOK && nTries < maxTries);
00317
00318
00319 if(inverseKinematics) {
00320 delete aProjectileTrack;
00321 delete theTargetNucleus;
00322 delete toInverseKinematics;
00323 delete toDirectKinematics;
00324 }
00325
00326 if(!eventIsOK) {
00327 std::stringstream ss;
00328 ss << "maximum number of tries exceeded for the proposed "
00329 << aTrack.GetKineticEnergy()/MeV << "-MeV " << aTrack.GetDefinition()->GetParticleName()
00330 << " + " << G4ParticleTable::GetParticleTable()->GetIon(theNucleus.GetZ_asInt(), theNucleus.GetA_asInt(), 0.0)->GetParticleName()
00331 << " inelastic reaction, in " << (inverseKinematics ? "inverse" : "direct") << " kinematics.";
00332 theInterfaceStore->EmitWarning(ss.str());
00333 theResult.SetStatusChange(isAlive);
00334 theResult.SetEnergyChange(aTrack.GetKineticEnergy());
00335 theResult.SetMomentumChange(aTrack.Get4Momentum().vect().unit());
00336 return &theResult;
00337 }
00338
00339
00340
00341 if(theExcitationHandler != 0) {
00342 for(std::list<G4Fragment>::const_iterator i = remnants.begin();
00343 i != remnants.end(); i++) {
00344 G4ReactionProductVector *deExcitationResult = theExcitationHandler->BreakItUp((*i));
00345
00346 for(G4ReactionProductVector::iterator fragment = deExcitationResult->begin();
00347 fragment != deExcitationResult->end(); ++fragment) {
00348 G4ParticleDefinition *def = (*fragment)->GetDefinition();
00349 if(def != 0) {
00350 G4DynamicParticle *theFragment = new G4DynamicParticle(def, (*fragment)->GetMomentum());
00351 theResult.AddSecondary(theFragment);
00352 }
00353 }
00354
00355 for(G4ReactionProductVector::iterator fragment = deExcitationResult->begin();
00356 fragment != deExcitationResult->end(); ++fragment) {
00357 delete (*fragment);
00358 }
00359 deExcitationResult->clear();
00360 delete deExcitationResult;
00361 }
00362 }
00363
00364 remnants.clear();
00365
00366 return &theResult;
00367 }
00368
00369 G4ReactionProductVector* G4INCLXXInterface::Propagate(G4KineticTrackVector* , G4V3DNucleus* ) {
00370 return 0;
00371 }
00372
00373 G4INCL::ParticleType G4INCLXXInterface::toINCLParticleType(G4ParticleDefinition const * const pdef) const {
00374 if( pdef == G4Proton::Proton()) return G4INCL::Proton;
00375 else if(pdef == G4Neutron::Neutron()) return G4INCL::Neutron;
00376 else if(pdef == G4PionPlus::PionPlus()) return G4INCL::PiPlus;
00377 else if(pdef == G4PionMinus::PionMinus()) return G4INCL::PiMinus;
00378 else if(pdef == G4PionZero::PionZero()) return G4INCL::PiZero;
00379 else if(pdef == G4Deuteron::Deuteron()) return G4INCL::Composite;
00380 else if(pdef == G4Triton::Triton()) return G4INCL::Composite;
00381 else if(pdef == G4He3::He3()) return G4INCL::Composite;
00382 else if(pdef == G4Alpha::Alpha()) return G4INCL::Composite;
00383 else if(pdef->GetParticleType() == G4GenericIon::GenericIon()->GetParticleType()) return G4INCL::Composite;
00384 else return G4INCL::UnknownParticle;
00385 }
00386
00387 G4INCL::ParticleSpecies G4INCLXXInterface::toINCLParticleSpecies(G4HadProjectile const &aTrack) const {
00388 const G4ParticleDefinition *pdef = aTrack.GetDefinition();
00389 const G4INCL::ParticleType theType = toINCLParticleType(pdef);
00390 if(theType!=G4INCL::Composite)
00391 return G4INCL::ParticleSpecies(theType);
00392 else {
00393 G4INCL::ParticleSpecies theSpecies;
00394 theSpecies.theType=theType;
00395 theSpecies.theA=(G4int) pdef->GetBaryonNumber();
00396 theSpecies.theZ=(G4int) pdef->GetPDGCharge();
00397 return theSpecies;
00398 }
00399 }
00400
00401 G4double G4INCLXXInterface::toINCLKineticEnergy(G4HadProjectile const &aTrack) const {
00402 return aTrack.GetKineticEnergy();
00403 }
00404
00405 G4ParticleDefinition *G4INCLXXInterface::toG4ParticleDefinition(G4int A, G4int Z) const {
00406 if (A == 1 && Z == 1) return G4Proton::Proton();
00407 else if(A == 1 && Z == 0) return G4Neutron::Neutron();
00408 else if(A == 0 && Z == 1) return G4PionPlus::PionPlus();
00409 else if(A == 0 && Z == -1) return G4PionMinus::PionMinus();
00410 else if(A == 0 && Z == 0) return G4PionZero::PionZero();
00411 else if(A == 2 && Z == 1) return G4Deuteron::Deuteron();
00412 else if(A == 3 && Z == 1) return G4Triton::Triton();
00413 else if(A == 3 && Z == 2) return G4He3::He3();
00414 else if(A == 4 && Z == 2) return G4Alpha::Alpha();
00415 else if(A > 0 && Z > 0 && A > Z) {
00416 return G4ParticleTable::GetParticleTable()->GetIon(Z, A, 0.0);
00417 } else {
00418 return 0;
00419 }
00420 }
00421
00422 G4DynamicParticle *G4INCLXXInterface::toG4Particle(G4int A, G4int Z,
00423 G4double kinE,
00424 G4double px,
00425 G4double py, G4double pz) const {
00426 const G4ParticleDefinition *def = toG4ParticleDefinition(A, Z);
00427 if(def == 0) {
00428 return 0;
00429 }
00430 const G4double energy = kinE / MeV;
00431 const G4ThreeVector momentum(px, py, pz);
00432 const G4ThreeVector momentumDirection = momentum.unit();
00433 G4DynamicParticle *p = new G4DynamicParticle(def, momentumDirection, energy);
00434 return p;
00435 }
00436
00437 G4double G4INCLXXInterface::remnant4MomentumScaling(G4double mass,
00438 G4double kineticE,
00439 G4double px, G4double py,
00440 G4double pz) const {
00441 const G4double p2 = px*px + py*py + pz*pz;
00442 if(p2 > 0.0) {
00443 const G4double pnew2 = kineticE*kineticE + 2.0*kineticE*mass;
00444 return std::sqrt(pnew2)/std::sqrt(p2);
00445 } else {
00446 return 1.0;
00447 }
00448 }
00449