G4QCoherentChargeExchange Class Reference

#include <G4QCoherentChargeExchange.hh>

Inheritance diagram for G4QCoherentChargeExchange:

G4VDiscreteProcess G4VProcess

Public Member Functions

 G4QCoherentChargeExchange (const G4String &processName="CHIPS_CoherChargeExScattering")
 ~G4QCoherentChargeExchange ()
G4bool IsApplicable (const G4ParticleDefinition &particle)
G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)
G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep)
G4LorentzVector GetEnegryMomentumConservation ()
G4int GetNumberOfNeutronsInTarget ()

Detailed Description

Definition at line 85 of file G4QCoherentChargeExchange.hh.


Constructor & Destructor Documentation

G4QCoherentChargeExchange::G4QCoherentChargeExchange ( const G4String processName = "CHIPS_CoherChargeExScattering"  ) 

Definition at line 75 of file G4QCoherentChargeExchange.cc.

References G4cout, G4endl, G4HadronicDeprecate, G4VProcess::GetProcessName(), and G4VProcess::verboseLevel.

00076   : G4VDiscreteProcess(processName, fHadronic)
00077 {
00078   G4HadronicDeprecate("G4QCoherentChargeExchange");
00079 
00080 #ifdef debug
00081   G4cout<<"G4QCohChargeEx::Constructor is called processName="<<processName<<G4endl;
00082 #endif
00083   if (verboseLevel>0) G4cout << GetProcessName() << " process is created "<< G4endl;
00084   //G4QCHIPSWorld::Get()->GetParticles(nPartCWorld); // Create CHIPS World (234 part. max)
00085 }

G4QCoherentChargeExchange::~G4QCoherentChargeExchange (  ) 

Definition at line 88 of file G4QCoherentChargeExchange.cc.

00088 {}


Member Function Documentation

G4LorentzVector G4QCoherentChargeExchange::GetEnegryMomentumConservation (  ) 

Definition at line 91 of file G4QCoherentChargeExchange.cc.

00092                                                                 {return EnMomConservation;}

G4double G4QCoherentChargeExchange::GetMeanFreePath ( const G4Track aTrack,
G4double  previousStepSize,
G4ForceCondition condition 
) [virtual]

Implements G4VDiscreteProcess.

Definition at line 99 of file G4QCoherentChargeExchange.cc.

References DBL_MAX, G4cerr, G4cout, G4endl, G4QIsotope::Get(), G4QIsotope::GetCSVector(), G4DynamicParticle::GetDefinition(), G4Track::GetDynamicParticle(), G4Material::GetElementVector(), G4DynamicParticle::GetKineticEnergy(), G4Track::GetMaterial(), G4QIsotope::GetMeanCrossSection(), G4Material::GetNumberOfElements(), G4DynamicParticle::GetTotalMomentum(), G4Material::GetVecNbOfAtomsPerVolume(), G4QIsotope::InitElement(), IsApplicable(), G4QIsotope::IsDefined(), G4Neutron::Neutron(), NotForced, and G4Proton::Proton().

Referenced by PostStepDoIt().

00101 {
00102   *Fc = NotForced;
00103   const G4DynamicParticle* incidentParticle = aTrack.GetDynamicParticle();
00104   G4ParticleDefinition* incidentParticleDefinition=incidentParticle->GetDefinition();
00105   if( !IsApplicable(*incidentParticleDefinition))
00106     G4cout<<"*W*G4QCohChargeEx::GetMeanFreePath called for notImplementedParticle"<<G4endl;
00107   // Calculate the mean Cross Section for the set of Elements(*Isotopes) in the Material
00108   G4double Momentum = incidentParticle->GetTotalMomentum(); // 3-momentum of the Particle
00109 #ifdef debug
00110   G4double KinEn = incidentParticle->GetKineticEnergy();
00111   G4cout<<"G4QCohChEx::GetMeanFreePath: kinE="<<KinEn<<",Mom="<<Momentum<<G4endl; // Result
00112 #endif
00113   const G4Material* material = aTrack.GetMaterial();        // Get the current material
00114   const G4double* NOfNucPerVolume = material->GetVecNbOfAtomsPerVolume();
00115   const G4ElementVector* theElementVector = material->GetElementVector();
00116   G4int nE=material->GetNumberOfElements();
00117 #ifdef debug
00118   G4cout<<"G4QCohChargeExchange::GetMeanFreePath:"<<nE<<" Elem's in theMaterial"<<G4endl;
00119 #endif
00120   G4int pPDG=0;
00121 
00122   if     (incidentParticleDefinition == G4Proton::Proton()  ) pPDG=2212;
00123   else if(incidentParticleDefinition == G4Neutron::Neutron()) pPDG=2112;
00124   else G4cout<<"G4QCohChargeEx::GetMeanFreePath: only nA & pA are implemented"<<G4endl;
00125   
00126   G4QIsotope* Isotopes = G4QIsotope::Get(); // Pointer to the G4QIsotopes singleton
00127   G4double sigma=0.;                        // Sums over elements for the material
00128   G4int IPIE=IsoProbInEl.size();            // How many old elements?
00129   if(IPIE) for(G4int ip=0; ip<IPIE; ++ip)   // Clean up the SumProb's of Isotopes (SPI)
00130   {
00131     std::vector<G4double>* SPI=IsoProbInEl[ip]; // Pointer to the SPI vector
00132     SPI->clear();
00133     delete SPI;
00134     std::vector<G4int>* IsN=ElIsoN[ip];     // Pointer to the N vector
00135     IsN->clear();
00136     delete IsN;
00137   }
00138   ElProbInMat.clear();                      // Clean up the SumProb's of Elements (SPE)
00139   ElementZ.clear();                         // Clear the body vector for Z of Elements
00140   IsoProbInEl.clear();                      // Clear the body vector for SPI
00141   ElIsoN.clear();                           // Clear the body vector for N of Isotopes
00142   for(G4int i=0; i<nE; ++i)
00143   {
00144     G4Element* pElement=(*theElementVector)[i]; // Pointer to the current element
00145     G4int Z = static_cast<G4int>(pElement->GetZ()); // Z of the Element
00146     ElementZ.push_back(Z);                  // Remember Z of the Element
00147     G4int isoSize=0;                        // The default for the isoVectorLength is 0
00148     G4int indEl=0;                          // Index of non-natural element or 0(default)
00149     G4IsotopeVector* isoVector=pElement->GetIsotopeVector(); // Get the predefined IsoVect
00150     if(isoVector) isoSize=isoVector->size();// Get size of the existing isotopeVector
00151 #ifdef debug
00152     G4cout<<"G4QCoherentChargeExchange::GetMeanFreePath:isovectorLength="<<isoSize<<G4endl;
00153 #endif
00154     if(isoSize)                             // The Element has non-trivial abundance set
00155     {
00156       indEl=pElement->GetIndex()+1;         // Index of the non-trivial element is an order
00157 #ifdef debug
00158       G4cout<<"G4QCCX::GetMFP: iE="<<indEl<<", def="<<Isotopes->IsDefined(Z,indEl)<<G4endl;
00159 #endif
00160       if(!Isotopes->IsDefined(Z,indEl))     // This index is not defined for this Z: define
00161       {
00162         std::vector<std::pair<G4int,G4double>*>* newAbund =
00163                                                new std::vector<std::pair<G4int,G4double>*>;
00164         G4double* abuVector=pElement->GetRelativeAbundanceVector();
00165         for(G4int j=0; j<isoSize; j++)      // Calculation of abundance vector for isotopes
00166         {
00167           G4int N=pElement->GetIsotope(j)->GetN()-Z; // N means A=N+Z !
00168           if(pElement->GetIsotope(j)->GetZ()!=Z)G4cerr<<"G4QCohChX::GetMeanFreePath: Z="
00169                                          <<pElement->GetIsotope(j)->GetZ()<<"#"<<Z<<G4endl;
00170           G4double abund=abuVector[j];
00171           std::pair<G4int,G4double>* pr= new std::pair<G4int,G4double>(N,abund);
00172 #ifdef debug
00173           G4cout<<"G4QCohChEx::GetMeanFreePath:pair#="<<j<<",N="<<N<<",ab="<<abund<<G4endl;
00174 #endif
00175           newAbund->push_back(pr);
00176         }
00177 #ifdef debug
00178         G4cout<<"G4QCohChEx::GetMeanFreePath: pairVectorLength="<<newAbund->size()<<G4endl;
00179 #endif
00180         indEl=G4QIsotope::Get()->InitElement(Z,indEl,newAbund); // definition of the newInd
00181         for(G4int k=0; k<isoSize; k++) delete (*newAbund)[k];   // Cleaning temporary
00182         delete newAbund; // Was "new" in the beginning of the name space
00183       }
00184     }
00185     std::vector<std::pair<G4int,G4double>*>* cs= Isotopes->GetCSVector(Z,indEl);//CSPointer
00186     std::vector<G4double>* SPI = new std::vector<G4double>; // Pointer to the SPI vector
00187     IsoProbInEl.push_back(SPI);
00188     std::vector<G4int>* IsN = new std::vector<G4int>; // Pointer to the N vector
00189     ElIsoN.push_back(IsN);
00190     G4int nIs=cs->size();                   // A#Of Isotopes in the Element
00191 #ifdef debug
00192     G4cout<<"G4QCohChargEx::GMFP:=***=>,#isot="<<nIs<<", Z="<<Z<<", indEl="<<indEl<<G4endl;
00193 #endif
00194     G4double susi=0.;                       // sum of CS over isotopes
00195     if(nIs) for(G4int j=0; j<nIs; j++)      // Calculate CS for eachIsotope of El
00196     {
00197       std::pair<G4int,G4double>* curIs=(*cs)[j]; // A pointer, which is used twice
00198       G4int N=curIs->first;                 // #of Neuterons in the isotope j of El i
00199       IsN->push_back(N);                    // Remember Min N for the Element
00200 #ifdef debug
00201       G4cout<<"G4QCCX::GMFP:true, P="<<Momentum<<",Z="<<Z<<",N="<<N<<",PDG="<<pPDG<<G4endl;
00202 #endif
00203       G4bool ccsf=true;
00204       if(Q==-27.) ccsf=false;
00205 #ifdef debug
00206       G4cout<<"G4QCoherentChargeExchange::GMFP: GetCS #1 j="<<j<<G4endl;
00207 #endif
00208       G4double CSI=CalculateXSt(ccsf, true, Momentum, Z, N, pPDG);// CS(j,i) for theIsotope
00209 
00210 #ifdef debug
00211       G4cout<<"G4QCohChEx::GetMeanFreePath:jI="<<j<<",Zt="<<Z<<",Nt="<<N<<",Mom="<<Momentum
00212             <<", XSec="<<CSI/millibarn<<G4endl;
00213 #endif
00214       curIs->second = CSI;
00215       susi+=CSI;                            // Make a sum per isotopes
00216       SPI->push_back(susi);                 // Remember summed cross-section
00217     } // End of temporary initialization of the cross sections in the G4QIsotope singeltone
00218     sigma+=Isotopes->GetMeanCrossSection(Z,indEl)*NOfNucPerVolume[i];//SUM(MeanCS*NOfNperV)
00219 #ifdef debug
00220     G4cout<<"G4QCohChEx::GMFP:<S>="<<Isotopes->GetMeanCrossSection(Z,indEl)<<",AddToSigma="
00221           <<Isotopes->GetMeanCrossSection(Z,indEl)*NOfNucPerVolume[i]<<G4endl;
00222 #endif
00223     ElProbInMat.push_back(sigma);
00224   } // End of LOOP over Elements
00225   // Check that cross section is not zero and return the mean free path
00226 #ifdef debug
00227   G4cout<<"G4QCoherentChargeExchange::GetMeanFreePath: MeanFreePath="<<1./sigma<<G4endl;
00228 #endif
00229   if(sigma > 0.) return 1./sigma;                 // Mean path [distance] 
00230   return DBL_MAX;
00231 }

G4int G4QCoherentChargeExchange::GetNumberOfNeutronsInTarget (  ) 

Definition at line 94 of file G4QCoherentChargeExchange.cc.

00094 {return nOfNeutrons;}

G4bool G4QCoherentChargeExchange::IsApplicable ( const G4ParticleDefinition particle  )  [virtual]

Reimplemented from G4VProcess.

Definition at line 233 of file G4QCoherentChargeExchange.cc.

References G4cout, G4endl, G4ParticleDefinition::GetPDGEncoding(), G4Neutron::Neutron(), and G4Proton::Proton().

Referenced by GetMeanFreePath(), and PostStepDoIt().

00234 {
00235   if      (particle == *(        G4Proton::Proton()        )) return true;
00236   else if (particle == *(       G4Neutron::Neutron()       )) return true;
00237   //else if (particle == *(     G4MuonMinus::MuonMinus()     )) return true; 
00238   //else if (particle == *(       G4TauPlus::TauPlus()       )) return true;
00239   //else if (particle == *(      G4TauMinus::TauMinus()      )) return true;
00240   //else if (particle == *(      G4Electron::Electron()      )) return true;
00241   //else if (particle == *(      G4Positron::Positron()      )) return true;
00242   //else if (particle == *(         G4Gamma::Gamma()         )) return true;
00243   //else if (particle == *(      G4MuonPlus::MuonPlus()      )) return true;
00244   //else if (particle == *(G4AntiNeutrinoMu::AntiNeutrinoMu())) return true;
00245   //else if (particle == *(   G4NeutrinoMu::NeutrinoMu()   )) return true;
00246   //else if (particle == *(     G4PionMinus::PionMinus()     )) return true;
00247   //else if (particle == *(      G4PionPlus::PionPlus()      )) return true;
00248   //else if (particle == *(      G4KaonPlus::KaonPlus()      )) return true;
00249   //else if (particle == *(     G4KaonMinus::KaonMinus()     )) return true;
00250   //else if (particle == *(  G4KaonZeroLong::KaonZeroLong()  )) return true;
00251   //else if (particle == *( G4KaonZeroShort::KaonZeroShort() )) return true;
00252   //else if (particle == *(        G4Lambda::Lambda()        )) return true;
00253   //else if (particle == *(     G4SigmaPlus::SigmaPlus()     )) return true;
00254   //else if (particle == *(    G4SigmaMinus::SigmaMinus()    )) return true;
00255   //else if (particle == *(     G4SigmaZero::SigmaZero()     )) return true;
00256   //else if (particle == *(       G4XiMinus::XiMinus()       )) return true;
00257   //else if (particle == *(        G4XiZero::XiZero()        )) return true;
00258   //else if (particle == *(    G4OmegaMinus::OmegaMinus()    )) return true;
00259   //else if (particle == *(   G4AntiNeutron::AntiNeutron()   )) return true;
00260   //else if (particle == *(    G4AntiProton::AntiProton()    )) return true;
00261 #ifdef debug
00262   G4cout<<"***>>G4QCoherChargeExch::IsApplicable: PDG="<<particle.GetPDGEncoding()<<G4endl;
00263 #endif
00264   return false;
00265 }

G4VParticleChange * G4QCoherentChargeExchange::PostStepDoIt ( const G4Track aTrack,
const G4Step aStep 
) [virtual]

Reimplemented from G4VDiscreteProcess.

Definition at line 267 of file G4QCoherentChargeExchange.cc.

References G4ParticleChange::AddSecondary(), G4VProcess::aParticleChange, G4ParticleTable::FindIon(), fStopAndKill, G4cerr, G4cout, G4endl, G4UniformRand, G4QCHIPSWorld::Get(), G4DynamicParticle::Get4Momentum(), G4DynamicParticle::GetDefinition(), G4Track::GetDynamicParticle(), G4Material::GetElementVector(), G4Track::GetGlobalTime(), G4DynamicParticle::GetKineticEnergy(), G4DynamicParticle::GetMass(), G4QPDGCode::GetMass(), G4Track::GetMaterial(), GetMeanFreePath(), G4DynamicParticle::GetMomentumDirection(), G4Material::GetNumberOfElements(), G4ParticleDefinition::GetParticleName(), G4QCHIPSWorld::GetParticles(), G4ParticleDefinition::GetParticleSubType(), G4ParticleTable::GetParticleTable(), G4ParticleDefinition::GetParticleType(), G4ParticleDefinition::GetPDGEncoding(), G4Track::GetPosition(), G4DynamicParticle::GetTotalMomentum(), G4Track::GetTouchableHandle(), G4Track::GetWeight(), G4Element::GetZ(), G4ParticleChange::Initialize(), IsApplicable(), G4Neutron::Neutron(), NotForced, position, G4VDiscreteProcess::PostStepDoIt(), G4ParticleChange::ProposeEnergy(), G4VParticleChange::ProposeLocalEnergyDeposit(), G4ParticleChange::ProposeMomentumDirection(), G4VParticleChange::ProposeTrackStatus(), G4Proton::Proton(), G4DynamicParticle::Set4Momentum(), G4DynamicParticle::SetDefinition(), G4Track::SetTouchableHandle(), and G4Track::SetWeight().

00269 {
00270   static const G4double mProt= G4QPDGCode(2212).GetMass();     // CHIPS pMass in MeV
00271   static const G4double mNeut= G4QPDGCode(2212).GetMass();     // CHIPS pMass in MeV
00272   //
00273   //-------------------------------------------------------------------------------------
00274   static G4bool CWinit = true;                       // CHIPS Warld needs to be initted
00275   if(CWinit)
00276   {
00277     CWinit=false;
00278     G4QCHIPSWorld::Get()->GetParticles(nPartCWorld); // Create CHIPS World (234 part.max)
00279   }
00280   //-------------------------------------------------------------------------------------
00281   const G4DynamicParticle* projHadron = track.GetDynamicParticle();
00282   const G4ParticleDefinition* particle=projHadron->GetDefinition();
00283 #ifdef debug
00284   G4cout<<"G4QCohChargeExchange::PostStepDoIt: Before the GetMeanFreePath is called In4M="
00285         <<projHadron->Get4Momentum()<<" of PDG="<<particle->GetPDGEncoding()<<", Type="
00286         <<particle->GetParticleType()<<", Subtp="<<particle->GetParticleSubType()<<G4endl;
00287 #endif
00288   G4ForceCondition cond=NotForced;
00289   GetMeanFreePath(track, -27., &cond);                  // @@ ?? jus to update parameters?
00290 #ifdef debug
00291   G4cout<<"G4QCohChargeExchange::PostStepDoIt: After GetMeanFreePath is called"<<G4endl;
00292 #endif
00293   G4LorentzVector proj4M=(projHadron->Get4Momentum())/MeV; // Convert to MeV!
00294   G4double momentum = projHadron->GetTotalMomentum()/MeV; // 3-momentum of the Proj in MeV
00295   G4double Momentum = proj4M.rho();                   // @@ Just for the test purposes
00296   if(std::fabs(Momentum-momentum)>.000001)
00297        G4cerr<<"*Warning*G4QCohChEx::PostStepDoIt:P(IU)="<<Momentum<<"#"<<momentum<<G4endl;
00298 #ifdef pdebug
00299   G4cout<<"G4QCoherentChargeExchange::PostStepDoIt: pP(IU)="<<Momentum<<"="<<momentum
00300         <<",pM="<<pM<<",proj4M="<<proj4M<<proj4M.m()<<G4endl;
00301 #endif
00302   if (!IsApplicable(*particle))  // Check applicability
00303   {
00304     G4cerr<<"G4QCoherentChargeExchange::PostStepDoIt: Only NA is implemented."<<G4endl;
00305     return 0;
00306   }
00307   const G4Material* material = track.GetMaterial();      // Get the current material
00308   G4int Z=0;
00309   const G4ElementVector* theElementVector = material->GetElementVector();
00310   G4int nE=material->GetNumberOfElements();
00311 #ifdef debug
00312   G4cout<<"G4QCohChargeExchange::PostStepDoIt: "<<nE<<" elements in the material."<<G4endl;
00313 #endif
00314   G4int projPDG=0;                           // PDG Code prototype for the captured hadron
00315   // Not all these particles are implemented yet (see Is Applicable)
00316   if      (particle ==          G4Proton::Proton()         ) projPDG= 2212;
00317   else if (particle ==         G4Neutron::Neutron()        ) projPDG= 2112;
00318   //else if (particle ==       G4PionMinus::PionMinus()      ) projPDG= -211;
00319   //else if (particle ==        G4PionPlus::PionPlus()       ) projPDG=  211;
00320   //else if (particle ==        G4KaonPlus::KaonPlus()       ) projPDG= 2112;
00321   //else if (particle ==       G4KaonMinus::KaonMinus()      ) projPDG= -321;
00322   //else if (particle ==    G4KaonZeroLong::KaonZeroLong()   ) projPDG=  130;
00323   //else if (particle ==   G4KaonZeroShort::KaonZeroShort()  ) projPDG=  310;
00324   //else if (particle ==        G4MuonPlus::MuonPlus()       ) projPDG=  -13;
00325   //else if (particle ==       G4MuonMinus::MuonMinus()      ) projPDG=   13;
00326   //else if (particle ==      G4NeutrinoMu::NeutrinoMu()     ) projPDG=   14;
00327   //else if (particle ==  G4AntiNeutrinoMu::AntiNeutrinoMu() ) projPDG=  -14;
00328   //else if (particle ==        G4Electron::Electron()       ) projPDG=   11;
00329   //else if (particle ==        G4Positron::Positron()       ) projPDG=  -11;
00330   //else if (particle ==       G4NeutrinoE::NeutrinoE()      ) projPDG=   12;
00331   //else if (particle ==   G4AntiNeutrinoE::AntiNeutrinoE()  ) projPDG=  -12;
00332   //else if (particle ==           G4Gamma::Gamma()          ) projPDG=   22;
00333   //else if (particle ==         G4TauPlus::TauPlus()        ) projPDG=  -15;
00334   //else if (particle ==        G4TauMinus::TauMinus()       ) projPDG=   15;
00335   //else if (particle ==     G4NeutrinoTau::NeutrinoTau()    ) projPDG=   16;
00336   //else if (particle == G4AntiNeutrinoTau::AntiNeutrinoTau()) projPDG=  -16;
00337   //else if (particle ==          G4Lambda::Lambda()         ) projPDG= 3122;
00338   //else if (particle ==       G4SigmaPlus::SigmaPlus()      ) projPDG= 3222;
00339   //else if (particle ==      G4SigmaMinus::SigmaMinus()     ) projPDG= 3112;
00340   //else if (particle ==       G4SigmaZero::SigmaZero()      ) projPDG= 3212;
00341   //else if (particle ==         G4XiMinus::XiMinus()        ) projPDG= 3312;
00342   //else if (particle ==          G4XiZero::XiZero()         ) projPDG= 3322;
00343   //else if (particle ==      G4OmegaMinus::OmegaMinus()     ) projPDG= 3334;
00344   //else if (particle ==     G4AntiNeutron::AntiNeutron()    ) projPDG=-2112;
00345   //else if (particle ==      G4AntiProton::AntiProton()     ) projPDG=-2212;
00346 #ifdef debug
00347   G4int prPDG=particle->GetPDGEncoding();
00348   G4cout<<"G4QCohChrgExchange::PostStepDoIt: projPDG="<<projPDG<<", stPDG="<<prPDG<<G4endl;
00349 #endif
00350   if(!projPDG)
00351   {
00352     G4cerr<<"*Warning*G4QCoherentChargeExchange::PostStepDoIt:UndefinedProjHadron"<<G4endl;
00353     return 0;
00354   }
00355   //G4double pM2=proj4M.m2();        // in MeV^2
00356   //G4double pM=std::sqrt(pM2);      // in MeV
00357   G4double pM=mNeut;
00358   if(projPDG==2112) pM=mProt;
00359   G4double pM2=pM*pM;
00360   // Element treatment
00361   G4int EPIM=ElProbInMat.size();
00362 #ifdef debug
00363   G4cout<<"G4QCohChEx::PostStDoIt:m="<<EPIM<<",n="<<nE<<",T="<<ElProbInMat[EPIM-1]<<G4endl;
00364 #endif
00365   G4int i=0;
00366   if(EPIM>1)
00367   {
00368     G4double rnd = ElProbInMat[EPIM-1]*G4UniformRand();
00369     for(i=0; i<nE; ++i)
00370     {
00371 #ifdef debug
00372       G4cout<<"G4QCohChEx::PostStepDoIt:EPM["<<i<<"]="<<ElProbInMat[i]<<",r="<<rnd<<G4endl;
00373 #endif
00374       if (rnd<ElProbInMat[i]) break;
00375     }
00376     if(i>=nE) i=nE-1;                        // Top limit for the Element
00377   }
00378   G4Element* pElement=(*theElementVector)[i];
00379   Z=static_cast<G4int>(pElement->GetZ());
00380 #ifdef debug
00381     G4cout<<"G4QCoherentChargeExchange::PostStepDoIt: i="<<i<<", Z(element)="<<Z<<G4endl;
00382 #endif
00383   if(Z<=0)
00384   {
00385     G4cerr<<"-Warning-G4QCoherentChargeExchange::PostStepDoIt: Element with Z="<<Z<<G4endl;
00386     if(Z<0) return 0;
00387   }
00388   std::vector<G4double>* SPI = IsoProbInEl[i];// Vector of summedProbabilities for isotopes
00389   std::vector<G4int>* IsN = ElIsoN[i];     // Vector of "#of neutrons" in the isotope El[i]
00390   G4int nofIsot=SPI->size();               // #of isotopes in the element i
00391 #ifdef debug
00392   G4cout<<"G4QCohChargeExchange::PosStDoIt:nI="<<nofIsot<<",T="<<(*SPI)[nofIsot-1]<<G4endl;
00393 #endif
00394   G4int j=0;
00395   if(nofIsot>1)
00396   {
00397     G4double rndI=(*SPI)[nofIsot-1]*G4UniformRand(); // Randomize the isotop of the Element
00398     for(j=0; j<nofIsot; ++j)
00399     {
00400 #ifdef debug
00401       G4cout<<"G4QCohChargEx::PostStepDoIt: SP["<<j<<"]="<<(*SPI)[j]<<", r="<<rndI<<G4endl;
00402 #endif
00403       if(rndI < (*SPI)[j]) break;
00404     }
00405     if(j>=nofIsot) j=nofIsot-1;            // Top limit for the isotope
00406   }
00407   G4int N =(*IsN)[j]; ;                    // Randomized number of neutrons
00408 #ifdef debug
00409   G4cout<<"G4QCohChargeEx::PostStepDoIt: j="<<i<<", N(isotope)="<<N<<", MeV="<<MeV<<G4endl;
00410 #endif
00411   if(N<0)
00412   {
00413     G4cerr<<"*Warning*G4QCohChEx::PostStepDoIt: Isotope with Z="<<Z<<", 0>N="<<N<<G4endl;
00414     return 0;
00415   }
00416   nOfNeutrons=N;                           // Remember it for the energy-momentum check
00417 #ifdef debug
00418   G4cout<<"G4QCohChargeExchange::PostStepDoIt: N="<<N<<" for element with Z="<<Z<<G4endl;
00419 #endif
00420   if(N<0)
00421   {
00422     G4cerr<<"*Warning*G4QCoherentChargeExchange::PostStepDoIt:Element with N="<<N<< G4endl;
00423     return 0;
00424   }
00425   aParticleChange.Initialize(track);
00426 #ifdef debug
00427   G4cout<<"G4QCoherentChargeExchange::PostStepDoIt: track is initialized"<<G4endl;
00428 #endif
00429   G4double      weight    = track.GetWeight();
00430   G4double      localtime = track.GetGlobalTime();
00431   G4ThreeVector position  = track.GetPosition();
00432 #ifdef debug
00433   G4cout<<"G4QCoherentChargeExchange::PostStepDoIt: before Touchable extraction"<<G4endl;
00434 #endif
00435   G4TouchableHandle trTouchable = track.GetTouchableHandle();
00436 #ifdef debug
00437   G4cout<<"G4QCoherentChargeExchange::PostStepDoIt: Touchable is extracted"<<G4endl;
00438 #endif
00439   //
00440   G4int targPDG=90000000+Z*1000+N;         // CHIPS PDG Code of the target nucleus
00441   if(projPDG==2212) targPDG+=999;          // convert to final nucleus code
00442   else if(projPDG==2112) targPDG-=999;
00443   G4QPDGCode targQPDG(targPDG);            // @@ use G4Ion and get rid of CHIPS World
00444   G4double tM=targQPDG.GetMass();          // CHIPS final nucleus mass in MeV
00445   G4double kinEnergy= projHadron->GetKineticEnergy()*MeV; // Kin energy in MeV (Is *MeV n?)
00446   G4ParticleMomentum dir = projHadron->GetMomentumDirection();// It is a unit three-vector
00447   G4LorentzVector tot4M=proj4M+G4LorentzVector(0.,0.,0.,tM); // Total 4-mom of the reaction
00448 #ifdef debug
00449   G4cout<<"G4QCohChrgEx::PostStepDoIt: tM="<<tM<<", p4M="<<proj4M<<", t4M="<<tot4M<<G4endl;
00450 #endif
00451   EnMomConservation=tot4M;                 // Total 4-mom of reaction for E/M conservation
00452   // @@ Probably this is not necessary any more
00453 #ifdef debug
00454   G4cout<<"G4QCCX::PSDI:false, P="<<Momentum<<",Z="<<Z<<",N="<<N<<",PDG="<<projPDG<<G4endl;
00455 #endif
00456   G4double xSec=CalculateXSt(false, true, Momentum, Z, N, projPDG); // Recalc. CrossSection
00457 #ifdef debug
00458   G4cout<<"G4QCoChEx::PSDI:PDG="<<projPDG<<",P="<<Momentum<<",CS="<<xSec/millibarn<<G4endl;
00459 #endif
00460 #ifdef nandebug
00461   if(xSec>0. || xSec<0. || xSec==0);
00462   else  G4cout<<"*Warning*G4QCohChargeExchange::PSDI: xSec="<<xSec/millibarn<<G4endl;
00463 #endif
00464   // @@ check a possibility to separate p, n, or alpha (!)
00465   if(xSec <= 0.) // The cross-section iz 0 -> Do Nothing
00466   {
00467 #ifdef pdebug
00468     G4cerr<<"*Warning*G4QCoherentChargeExchange::PSDoIt:*Zero cross-section* PDG="<<projPDG
00469           <<",tPDG="<<targPDG<<",P="<<Momentum<<G4endl;
00470 #endif
00471     //Do Nothing Action insead of the reaction
00472     aParticleChange.ProposeEnergy(kinEnergy);
00473     aParticleChange.ProposeLocalEnergyDeposit(0.);
00474     aParticleChange.ProposeMomentumDirection(dir) ;
00475     return G4VDiscreteProcess::PostStepDoIt(track,step);
00476   }
00477   G4double mint=CalculateXSt(false, false, Momentum, Z, N, projPDG);// randomize t in MeV^2
00478 #ifdef pdebug
00479   G4cout<<"G4QCohChEx::PoStDoIt:pPDG="<<projPDG<<",tPDG="<<targPDG<<",P="<<Momentum<<",CS="
00480         <<xSec<<",-t="<<mint<<G4endl;
00481 #endif
00482 #ifdef nandebug
00483   if(mint>-.0000001);
00484   else  G4cout<<"*Warning*G4QCoherentChargeExchange::PostStDoIt:-t="<<mint<<G4endl;
00485 #endif
00486   G4double maxt=CalculateXSt(true, false, Momentum, Z, N, projPDG);
00487   if(maxt<=0.) maxt=1.e-27;
00488   G4double cost=1.-mint/maxt; // cos(theta) in CMS (@@ we neglect mass diff for ChEx)
00489   // 
00490 #ifdef ppdebug
00491   G4cout<<"G4QCoherentChargeExchange::PoStDoIt:t="<<mint<<",dpcm2="<<maxt
00492         <<",Ek="<<kinEnergy<<",tM="<<tM<<",pM="<<pM<<",cost="<<cost<<G4endl;
00493 #endif
00494   if(cost>1. || cost<-1. || !(cost>-1. || cost<1.))
00495   {
00496     if(cost>1.000001 || cost<-1.000001 || !(cost>-1. || cost<1.))
00497     {
00498       G4double tM2=tM*tM;                         // Squared target mass
00499       G4double pEn=pM+kinEnergy;                  // tot projectile Energy in MeV
00500       G4double sM=(tM+tM)*pEn+tM2+pM2;            // Mondelstam s
00501       G4double twop2cm=(tM2+tM2)*(pEn*pEn-pM2)/sM;// Max_t/2 (2*p^2_cm)
00502       G4cout<<"*Warning*G4QCohChEx::PostStepDoIt:cos="<<cost<<",t="<<mint<<",T="<<kinEnergy
00503             <<",tM="<<tM<<",tmax="<<2*kinEnergy*tM<<",p="<<projPDG<<",t="<<targPDG<<G4endl;
00504       G4cout<<"..G4QCohChEx::PoStDoI: dpcm2="<<twop2cm<<"="<<maxt<<G4endl;
00505     }
00506     if     (cost>1.)  cost=1.;
00507     else if(cost<-1.) cost=-1.;
00508   }
00509   G4LorentzVector reco4M=G4LorentzVector(0.,0.,0.,tM);      // 4mom of the recoil target
00510   G4LorentzVector scat4M=G4LorentzVector(0.,0.,0.,pM);      // 4mom of the recoil target
00511   G4LorentzVector dir4M=tot4M-G4LorentzVector(0.,0.,0.,(tot4M.e()-tM-pM)*.01);
00512   if(!G4QHadron(tot4M).RelDecayIn2(scat4M, reco4M, dir4M, cost, cost))
00513   {
00514     G4cerr<<"G4QCohChEx::PSDI:t4M="<<tot4M<<",pM="<<pM<<",tM="<<tM<<",cost="<<cost<<G4endl;
00515     //G4Exception("G4QCoherentChargeExchange::PostStepDoIt:","009",FatalException,"Decay");
00516   }
00517 #ifdef debug
00518   G4cout<<"G4QCohChEx::PoStDoIt:s4M="<<scat4M<<"+r4M="<<reco4M<<"="<<scat4M+reco4M<<G4endl;
00519   G4cout<<"G4QCohChEx::PoStDoIt: scatE="<<scat4M.e()-pM<<", recoE="<<reco4M.e()-tM<<",d4M="
00520         <<tot4M-scat4M-reco4M<<G4endl;
00521 #endif
00522   // Kill scattered hadron
00523   aParticleChange.ProposeTrackStatus(fStopAndKill);
00524   // Definition of the scattered nucleon
00525   G4DynamicParticle* theSec = new G4DynamicParticle; // A secondary for the recoil hadron 
00526   G4ParticleDefinition* theDefinition=G4Proton::Proton();
00527   if(projPDG==2212) theDefinition=G4Neutron::Neutron();
00528   theSec->SetDefinition(theDefinition);
00529   EnMomConservation-=scat4M;
00530   theSec->Set4Momentum(scat4M);
00531   G4Track* aNewTrack = new G4Track(theSec, localtime, position );
00532   aNewTrack->SetWeight(weight);                                   //    weighted
00533   aNewTrack->SetTouchableHandle(trTouchable);
00534   aParticleChange.AddSecondary( aNewTrack );
00535   // Filling the recoil nucleus
00536   theSec = new G4DynamicParticle; // A secondary for the recoil hadron 
00537   G4int aA = Z+N;
00538 #ifdef pdebug
00539   G4cout<<"G4QCoherentChargeExchange::PostStepDoIt: Ion Z="<<Z<<", A="<<aA<<G4endl;
00540 #endif
00541   theDefinition=G4ParticleTable::GetParticleTable()->FindIon(Z,aA,0,Z);
00542   if(!theDefinition)G4cout<<"*Warning*G4QCohChEx::PostStepDoIt:drop PDG="<<targPDG<<G4endl;
00543 #ifdef pdebug
00544   G4cout<<"G4QCohChEx::PostStepDoIt:RecoilName="<<theDefinition->GetParticleName()<<G4endl;
00545 #endif
00546   theSec->SetDefinition(theDefinition);
00547   EnMomConservation-=reco4M;
00548 #ifdef tdebug
00549   G4cout<<"G4QCohChEx::PostSDoIt:"<<targPDG<<reco4M<<reco4M.m()<<EnMomConservation<<G4endl;
00550 #endif
00551   theSec->Set4Momentum(reco4M);
00552 #ifdef debug
00553   G4ThreeVector curD=theSec->GetMomentumDirection();
00554   G4double curM=theSec->GetMass();
00555   G4double curE=theSec->GetKineticEnergy()+curM;
00556   G4cout<<"G4QCohChEx::PostStpDoIt:p="<<curD<<curD.mag()<<",e="<<curE<<",m="<<curM<<G4endl;
00557 #endif
00558   // Make a recoil nucleus
00559   aNewTrack = new G4Track(theSec, localtime, position );
00560   aNewTrack->SetWeight(weight);                                   //    weighted
00561   aNewTrack->SetTouchableHandle(trTouchable);
00562   aParticleChange.AddSecondary( aNewTrack );
00563 #ifdef debug
00564     G4cout<<"G4QCoherentChargeExchange::PostStepDoIt:*** PostStepDoIt is done ***"<<G4endl;
00565 #endif
00566   return G4VDiscreteProcess::PostStepDoIt(track, step);
00567 }


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:53:05 2013 for Geant4 by  doxygen 1.4.7