// // ******************************************************************** // * License and Disclaimer * // * * // * The Geant4 software is copyright of the Copyright Holders of * // * the Geant4 Collaboration. It is provided under the terms and * // * conditions of the Geant4 Software License, included in the file * // * LICENSE and available at http://cern.ch/geant4/license . These * // * include a list of copyright holders. * // * * // * Neither the authors of this software system, nor their employing * // * institutes,nor the agencies providing financial support for this * // * work make any representation or warranty, express or implied, * // * regarding this software system or assume any liability for its * // * use. Please see the license in the file LICENSE and URL above * // * for the full disclaimer and the limitation of liability. * // * * // * This code implementation is the result of the scientific and * // * technical work of the GEANT4 collaboration. * // * By using, copying, modifying or distributing the software (or * // * any work based on the software) you agree to acknowledge its * // * use in resulting scientific publications, and indicate your * // * acceptance of all terms of the Geant4 Software license. * // ******************************************************************** // // $Id: G4QIonIonElastic.cc,v 1.4 2009/02/23 09:49:24 mkossov Exp $ // GEANT4 tag $Name: geant4-09-03-beta-cand-01 $ // // ---------------- G4QIonIonElastic class ----------------- // by Mikhail Kossov, December 2006. // G4QIonIonElastic class of the CHIPS Simulation Branch in GEANT4 // --------------------------------------------------------------- // **************************************************************************************** // ********** This CLASS is temporary moved from the photolepton_hadron directory ********* // **************************************************************************************** // Short description: a simple process for the Ion-Ion elastic scattering. // For heavy by heavy ions it can reach 50% of the total cross-section. // ----------------------------------------------------------------------- //#define debug //#define pdebug //#define tdebug //#define nandebug //#define ppdebug #include "G4QIonIonElastic.hh" // Initialization of static vectors G4int G4QIonIonElastic::nPartCWorld=152; // The#of particles initialized in CHIPS World std::vector G4QIonIonElastic::ElementZ; // Z of the element(i) in theLastCalc std::vector G4QIonIonElastic::ElProbInMat; // SumProbabilityElements in Material std::vector*> G4QIonIonElastic::ElIsoN; // N of isotope(j) of Element(i) std::vector*>G4QIonIonElastic::IsoProbInEl;//SumProbabIsotopes in ElI // Constructor G4QIonIonElastic::G4QIonIonElastic(const G4String& processName): G4VDiscreteProcess(processName, fHadronic) { #ifdef debug G4cout<<"G4QIonIonElastic::Constructor is called processName="<0) G4cout << GetProcessName() << " process is created "<< G4endl; SetProcessSubType(fHadronElastic); //G4QCHIPSWorld::Get()->GetParticles(nPartCWorld); // Create CHIPS World (234 part. max) } // Destructor G4QIonIonElastic::~G4QIonIonElastic() {} // output of the function must be in units of length! L=1/sig_V,sig_V=SUM(n(j,i)*sig(j,i)), // where n(i,j) is a number of nuclei of the isotop j of the element i in V=1(lengtUnit^3) // ********** All CHIPS cross sections are calculated in the surface units ************ G4double G4QIonIonElastic::GetMeanFreePath(const G4Track& aTrack, G4double, G4ForceCondition* Fc) { *Fc = NotForced; const G4DynamicParticle* incidentParticle = aTrack.GetDynamicParticle(); G4ParticleDefinition* incidentParticleDefinition=incidentParticle->GetDefinition(); if( !IsApplicable(*incidentParticleDefinition)) G4cout<<"-Warning-G4QIonIonElastic::GetMeanFreePath: notImplementedParticle"<GetTotalMomentum(); // 3-momentum of the Particle #ifdef debug G4double KinEn = incidentParticle->GetKineticEnergy(); G4cout<<"G4QIonIonElastic::GetMeanFreePath: kinE="<GetVecNbOfAtomsPerVolume(); const G4ElementVector* theElementVector = material->GetElementVector(); G4int nE=material->GetNumberOfElements(); #ifdef debug G4cout<<"G4QIonIonElastic::GetMeanFreePath:"<GetPDGEncoding(); if ( incidentParticleDefinition == G4Deuteron::Deuteron() ) pPDG = 100001002; else if ( incidentParticleDefinition == G4Alpha::Alpha() ) pPDG = 100002004; else if ( incidentParticleDefinition == G4Triton::Triton() ) pPDG = 100001003; else if ( incidentParticleDefinition == G4He3::He3() ) pPDG = 100002003; else if ( incidentParticleDefinition == G4GenericIon::GenericIon() ) { pPDG=incidentParticleDefinition->GetPDGEncoding(); #ifdef debug G4int B=incidentParticleDefinition->GetBaryonNumber(); G4int C=incidentParticleDefinition->GetPDGCharge(); prPDG=100000000+1000*C+B; G4cout<<"G4QIonIonElastic::GetMeanFreePath: PDG="<GetBaryonNumber(); // Divide Mom by projectile A G4QIsotope* Isotopes = G4QIsotope::Get(); // Pointer to the G4QIsotopes singleton G4double sigma=0.; // Sums over elements for the material G4int IPIE=IsoProbInEl.size(); // How many old elements? if(IPIE) for(G4int ip=0; ip* SPI=IsoProbInEl[ip]; // Pointer to the SPI vector SPI->clear(); delete SPI; std::vector* IsN=ElIsoN[ip]; // Pointer to the N vector IsN->clear(); delete IsN; } ElProbInMat.clear(); // Clean up the SumProb's of Elements (SPE) ElementZ.clear(); // Clear the body vector for Z of Elements IsoProbInEl.clear(); // Clear the body vector for SPI ElIsoN.clear(); // Clear the body vector for N of Isotopes for(G4int i=0; i(pElement->GetZ()); // Z of the Element ElementZ.push_back(Z); // Remember Z of the Element G4int isoSize=0; // The default for the isoVectorLength is 0 G4int indEl=0; // Index of non-natural element or 0(default) G4IsotopeVector* isoVector=pElement->GetIsotopeVector(); // Get the predefined IsoVect if(isoVector) isoSize=isoVector->size();// Get size of the existing isotopeVector #ifdef debug G4cout<<"G4QIonIonElastic::GetMeanFreePath: isovectorLength="<GetIndex()+1; // Index of the non-trivial element is an order #ifdef debug G4cout<<"G4QIIEl::GetMFP:iE="<* pr= new std::pair(N,abund); #ifdef debug G4cout<<"G4QIonIonElastic::GetMeanFP:pair#="<push_back(pr); } #ifdef debug G4cout<<"G4QIonIonElastic::GetMeanFP: pairVectorLength="<size()<InitElement(Z,indEl,newAbund); // definition of the newInd for(G4int k=0; k*>* cs= Isotopes->GetCSVector(Z,indEl);//CSPointer std::vector* SPI = new std::vector; // Pointer to the SPI vector IsoProbInEl.push_back(SPI); std::vector* IsN = new std::vector; // Pointer to the N vector ElIsoN.push_back(IsN); G4int nIs=cs->size(); // A#Of Isotopes in the Element #ifdef debug G4cout<<"G4QIonIonEl::GetMFP:=***=> #isot="<="<GetMeanCrossSection(Z,indEl)<<",AddToSig=" <GetMeanCrossSection(Z,indEl)*NOfNucPerVolume[i]<>G4QIonIonElastic::IsApplicable: PDG="<GetParticles(nPartCWorld); // Create CHIPS World (234 part.max) } //------------------------------------------------------------------------------------- const G4DynamicParticle* projHadron = track.GetDynamicParticle(); const G4ParticleDefinition* particle=projHadron->GetDefinition(); #ifdef debug G4cout<<"G4QIonIonElastic::PostStepDoIt: Before the GetMeanFreePath is called In4M=" <Get4Momentum()<<" of PDG="<GetPDGEncoding()<<", Type=" <GetParticleType()<<", Subtp="<GetParticleSubType()<Get4Momentum())/MeV; // Convert to MeV! G4LorentzVector scat4M=proj4M; // @@ Must be filled (?) G4double momentum = projHadron->GetTotalMomentum()/MeV; // 3-momentum of the Proj in MeV G4double Momentum = proj4M.rho(); // @@ Just for the test purposes if(std::fabs(Momentum-momentum)>.000001) G4cerr<<"-Warn-G4QIonIonElastic::PostStepDoIt:P(IU)="<GetElementVector(); G4int nE=material->GetNumberOfElements(); #ifdef debug G4cout<<"G4QIonIonElastic::PostStepDoIt: "<GetPDGEncoding(); G4int projPDG=0; // CHIPS PDG Code for the captured hadron if (particle == G4Deuteron::Deuteron() ) projPDG= 100001002; else if (particle == G4Alpha::Alpha() ) projPDG= 100002004; else if (particle == G4Triton::Triton() ) projPDG= 100001003; else if (particle == G4He3::He3() ) projPDG= 100002003; else if (particle == G4GenericIon::GenericIon() ) { projPDG=particle->GetPDGEncoding(); #ifdef debug G4int B=particle->GetBaryonNumber(); G4int C=particle->GetPDGCharge(); prPDG=100000000+1000*C+B; G4cout<<"G4QIonIonElastic::PostStepDoIt: PDG="<GetPDGEncoding(); G4cout<<"G4QIonIonElastic::PostStepDoIt: projPDG="<GetBaryonNumber(); // Projectile A G4double pZ=particle->GetPDGCharge(); // Projectile Z G4double pN=pA-pZ; // Projectile N G4int EPIM=ElProbInMat.size(); #ifdef debug G4cout<<"G4QIonIonElastic::PSDI:m="<1) { G4double rnd = ElProbInMat[EPIM-1]*G4UniformRand(); for(i=0; isize(); // #of isotopes in the element i #ifdef debug G4cout<<"G4QIonIonElastic::PosStDoIt: nI="<1) { G4double rndI=(*SPI)[nofIsot-1]*G4UniformRand(); // Randomize the isotop of the Element for(j=0; jN="<GetKineticEnergy()*MeV; // Kin energy in MeV (Is *MeV n?) G4ParticleMomentum dir = projHadron->GetMomentumDirection();// It is a unit three-vector G4LorentzVector tot4M=proj4M+G4LorentzVector(0.,0.,0.,tM); // Total 4-mom of the reaction #ifdef debug G4cout<<"G4QIonIonElastic::PostStDI: tM="<GetHMaxT()<1.) cost=1.; else if(cost<-1.) cost=-1.; } G4LorentzVector reco4M=G4LorentzVector(0.,0.,0.,tM); // 4mom of the recoil target G4LorentzVector dir4M=tot4M-G4LorentzVector(0.,0.,0.,(tot4M.e()-tM-pM)*.01); if(!G4QHadron(tot4M).RelDecayIn2(scat4M, reco4M, dir4M, cost, cost)) { G4cerr<<"G4QIonIonE::PSDI:t4M="<Set4Momentum(reco4M); #ifdef debug G4ThreeVector curD=theSec->GetMomentumDirection(); G4double curM=theSec->GetMass(); G4double curE=theSec->GetKineticEnergy()+curM; G4cout<<"G4QIonIonElastic::PSDI: p="<SetWeight(weight); // weighted aNewTrack->SetTouchableHandle(trTouchable); aParticleChange.AddSecondary( aNewTrack ); #ifdef debug G4cout<<"G4QIonIonElastic::PostStepDoIt: **** PostStepDoIt is done ****"<