G4RPGPiPlusInelastic.cc

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 // $Id$
00027 //
00028  
00029 #include "G4RPGPiPlusInelastic.hh"
00030 #include "G4SystemOfUnits.hh"
00031 #include "Randomize.hh"
00032 
00033 G4HadFinalState*
00034 G4RPGPiPlusInelastic::ApplyYourself(const G4HadProjectile& aTrack,
00035                                      G4Nucleus& targetNucleus)
00036 {
00037   const G4HadProjectile *originalIncident = &aTrack;
00038   if (originalIncident->GetKineticEnergy()<= 0.1) {
00039     theParticleChange.SetStatusChange(isAlive);
00040     theParticleChange.SetEnergyChange(aTrack.GetKineticEnergy());
00041     theParticleChange.SetMomentumChange(aTrack.Get4Momentum().vect().unit()); 
00042     return &theParticleChange;      
00043   }
00044 
00045     // create the target particle
00046     
00047     G4DynamicParticle *originalTarget = targetNucleus.ReturnTargetParticle();
00048     G4ReactionProduct targetParticle( originalTarget->GetDefinition() );
00049     
00050     G4ReactionProduct currentParticle( 
00051     const_cast<G4ParticleDefinition *>(originalIncident->GetDefinition() ) );
00052     currentParticle.SetMomentum( originalIncident->Get4Momentum().vect() );
00053     currentParticle.SetKineticEnergy( originalIncident->GetKineticEnergy() );
00054     
00055     // Fermi motion and evaporation
00056     // As of Geant3, the Fermi energy calculation had not been Done
00057     
00058     G4double ek = originalIncident->GetKineticEnergy();
00059     G4double amas = originalIncident->GetDefinition()->GetPDGMass();
00060     
00061     G4double tkin = targetNucleus.Cinema( ek );
00062     ek += tkin;
00063     currentParticle.SetKineticEnergy( ek );
00064     G4double et = ek + amas;
00065     G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) );
00066     G4double pp = currentParticle.GetMomentum().mag();
00067     if( pp > 0.0 )
00068     {
00069       G4ThreeVector momentum = currentParticle.GetMomentum();
00070       currentParticle.SetMomentum( momentum * (p/pp) );
00071     }
00072     
00073     // calculate black track energies
00074     
00075     tkin = targetNucleus.EvaporationEffects( ek );
00076     ek -= tkin;
00077     currentParticle.SetKineticEnergy( ek );
00078     et = ek + amas;
00079     p = std::sqrt( std::abs((et-amas)*(et+amas)) );
00080     pp = currentParticle.GetMomentum().mag();
00081     if( pp > 0.0 )
00082     {
00083       G4ThreeVector momentum = currentParticle.GetMomentum();
00084       currentParticle.SetMomentum( momentum * (p/pp) );
00085     }
00086 
00087     G4ReactionProduct modifiedOriginal = currentParticle;
00088 
00089     currentParticle.SetSide( 1 ); // incident always goes in forward hemisphere
00090     targetParticle.SetSide( -1 );  // target always goes in backward hemisphere
00091     G4bool incidentHasChanged = false;
00092     G4bool targetHasChanged = false;
00093     G4bool quasiElastic = false;
00094     G4FastVector<G4ReactionProduct,256> vec;  // vec will contain the secondary particles
00095     G4int vecLen = 0;
00096     vec.Initialize( 0 );
00097     
00098     const G4double cutOff = 0.1;
00099     if( currentParticle.GetKineticEnergy() > cutOff )
00100       InitialCollision(vec, vecLen, currentParticle, targetParticle,
00101                        incidentHasChanged, targetHasChanged);
00102     
00103     CalculateMomenta( vec, vecLen,
00104                       originalIncident, originalTarget, modifiedOriginal,
00105                       targetNucleus, currentParticle, targetParticle,
00106                       incidentHasChanged, targetHasChanged, quasiElastic );
00107     
00108     SetUpChange( vec, vecLen,
00109                  currentParticle, targetParticle,
00110                  incidentHasChanged );
00111     
00112     delete originalTarget;
00113     return &theParticleChange;
00114 }
00115 
00116 
00117 // Initial Collision
00118 //   selects the particle types arising from the initial collision of
00119 //   the projectile and target nucleon.  Secondaries are assigned to 
00120 //   forward and backward reaction hemispheres, but final state energies
00121 //   and momenta are not calculated here.
00122  
00123 void 
00124 G4RPGPiPlusInelastic::InitialCollision(G4FastVector<G4ReactionProduct,256>& vec,
00125                                   G4int& vecLen,
00126                                   G4ReactionProduct& currentParticle,
00127                                   G4ReactionProduct& targetParticle,
00128                                   G4bool& incidentHasChanged,
00129                                   G4bool& targetHasChanged)
00130 {
00131   G4double KE = currentParticle.GetKineticEnergy()/GeV;
00132 
00133   G4int mult;
00134   G4int partType;
00135   std::vector<G4int> fsTypes;
00136 
00137   G4double testCharge;
00138   G4double testBaryon;
00139   G4double testStrange;
00140 
00141   // Get particle types according to incident and target types
00142 
00143   if (targetParticle.GetDefinition() == particleDef[pro]) {
00144     mult = GetMultiplicityT32(KE);
00145     fsTypes = GetFSPartTypesForPipP(mult, KE);
00146     partType = fsTypes[0];
00147     if (partType != pro) {
00148       targetHasChanged = true;
00149       targetParticle.SetDefinition(particleDef[partType]);
00150     }
00151 
00152     testCharge = 2.0;
00153     testBaryon = 1.0;
00154     testStrange = 0.0;
00155 
00156   } else {   // target was a neutron
00157     mult = GetMultiplicityT12(KE);
00158     fsTypes = GetFSPartTypesForPipN(mult, KE);
00159     partType = fsTypes[0];
00160     if (partType != neu) {
00161       targetHasChanged = true;
00162       targetParticle.SetDefinition(particleDef[partType]);
00163     }
00164 
00165     testCharge = 1.0;
00166     testBaryon = 1.0;
00167     testStrange = 0.0;
00168   }
00169 
00170   // Remove target particle from list
00171 
00172   fsTypes.erase(fsTypes.begin());
00173 
00174   // See if the incident particle changed type 
00175 
00176   G4int choose = -1;
00177   for(G4int i=0; i < mult-1; ++i ) {
00178     partType = fsTypes[i];
00179     if (partType == pip) {
00180       choose = i;
00181       break;
00182     }
00183   }
00184   if (choose == -1) {
00185     incidentHasChanged = true;
00186     choose = G4int(G4UniformRand()*(mult-1) );
00187     partType = fsTypes[choose];
00188     currentParticle.SetDefinition(particleDef[partType]);
00189   }
00190   fsTypes.erase(fsTypes.begin()+choose);
00191 
00192   // Remaining particles are secondaries.  Put them into vec.
00193   //   Improve this by randomizing secondary order, then alternate
00194   //   which secondary is put into forward or backward hemisphere
00195 
00196   G4ReactionProduct* rp(0);
00197   for(G4int i=0; i < mult-2; ++i ) {
00198     partType = fsTypes[i];
00199     rp = new G4ReactionProduct();
00200     rp->SetDefinition(particleDef[partType]);
00201     (G4UniformRand() < 0.5) ? rp->SetSide(-1) : rp->SetSide(1);
00202     if (partType > pim && partType < pro) rp->SetMayBeKilled(false);  // kaons
00203     vec.SetElement(vecLen++, rp);
00204   }
00205  
00206   //  if (mult == 2 && !incidentHasChanged && !targetHasChanged) 
00207   //                                              quasiElastic = true;
00208 
00209   // Check conservation of charge, strangeness, baryon number
00210 
00211   CheckQnums(vec, vecLen, currentParticle, targetParticle,
00212              testCharge, testBaryon, testStrange);
00213 
00214   return;
00215 }

Generated on Mon May 27 17:49:45 2013 for Geant4 by  doxygen 1.4.7