source: trunk/source/processes/hadronic/models/chiral_inv_phase_space/interface/src/G4QIonIonElastic.cc @ 1007

Last change on this file since 1007 was 1007, checked in by garnier, 15 years ago

update to geant4.9.2

File size: 24.9 KB
Line 
1//
2// ********************************************************************
3// * License and Disclaimer                                           *
4// *                                                                  *
5// * The  Geant4 software  is  copyright of the Copyright Holders  of *
6// * the Geant4 Collaboration.  It is provided  under  the terms  and *
7// * conditions of the Geant4 Software License,  included in the file *
8// * LICENSE and available at  http://cern.ch/geant4/license .  These *
9// * include a list of copyright holders.                             *
10// *                                                                  *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work  make  any representation or  warranty, express or implied, *
14// * regarding  this  software system or assume any liability for its *
15// * use.  Please see the license in the file  LICENSE  and URL above *
16// * for the full disclaimer and the limitation of liability.         *
17// *                                                                  *
18// * This  code  implementation is the result of  the  scientific and *
19// * technical work of the GEANT4 collaboration.                      *
20// * By using,  copying,  modifying or  distributing the software (or *
21// * any work based  on the software)  you  agree  to acknowledge its *
22// * use  in  resulting  scientific  publications,  and indicate your *
23// * acceptance of all terms of the Geant4 Software license.          *
24// ********************************************************************
25//
26// $Id: G4QIonIonElastic.cc,v 1.3 2008/10/02 21:10:07 dennis Exp $
27// GEANT4 tag $Name: geant4-09-02 $
28//
29//      ---------------- G4QIonIonElastic class -----------------
30//                 by Mikhail Kossov, December 2006.
31// G4QIonIonElastic class of the CHIPS Simulation Branch in GEANT4
32// ---------------------------------------------------------------
33// ****************************************************************************************
34// ********** This CLASS is temporary moved from the photolepton_hadron directory *********
35// ****************************************************************************************
36
37//#define debug
38//#define pdebug
39//#define tdebug
40//#define nandebug
41//#define ppdebug
42
43#include "G4QIonIonElastic.hh"
44
45// Initialization of static vectors
46G4int G4QIonIonElastic::nPartCWorld=152;     // The#of particles initialized in CHIPS World
47std::vector<G4int> G4QIonIonElastic::ElementZ;        // Z of the element(i) in theLastCalc
48std::vector<G4double> G4QIonIonElastic::ElProbInMat;  // SumProbabilityElements in Material
49std::vector<std::vector<G4int>*> G4QIonIonElastic::ElIsoN; // N of isotope(j) of Element(i)
50std::vector<std::vector<G4double>*>G4QIonIonElastic::IsoProbInEl;//SumProbabIsotopes in ElI
51
52// Constructor
53G4QIonIonElastic::G4QIonIonElastic(const G4String& processName):
54  G4VDiscreteProcess(processName, fHadronic)
55{
56#ifdef debug
57  G4cout<<"G4QIonIonElastic::Constructor is called processName="<<processName<<G4endl;
58#endif
59  if (verboseLevel>0) G4cout << GetProcessName() << " process is created "<< G4endl;
60  SetProcessSubType(fHadronElastic);
61  //G4QCHIPSWorld::Get()->GetParticles(nPartCWorld); // Create CHIPS World (234 part. max)
62}
63
64// Destructor
65G4QIonIonElastic::~G4QIonIonElastic() {}
66
67// output of the function must be in units of length! L=1/sig_V,sig_V=SUM(n(j,i)*sig(j,i)),
68// where n(i,j) is a number of nuclei of the isotop j of the element i in V=1(lengtUnit^3)
69// ********** All CHIPS cross sections are calculated in the surface units ************
70G4double G4QIonIonElastic::GetMeanFreePath(const G4Track& aTrack, G4double,
71                                           G4ForceCondition* Fc)
72{
73  *Fc = NotForced;
74  const G4DynamicParticle* incidentParticle = aTrack.GetDynamicParticle();
75  G4ParticleDefinition* incidentParticleDefinition=incidentParticle->GetDefinition();
76  if( !IsApplicable(*incidentParticleDefinition))
77    G4cout<<"-Warning-G4QIonIonElastic::GetMeanFreePath: notImplementedParticle"<<G4endl;
78  // Calculate the mean Cross Section for the set of Elements(*Isotopes) in the Material
79  G4double Momentum = incidentParticle->GetTotalMomentum(); // 3-momentum of the Particle
80#ifdef debug
81  G4double KinEn = incidentParticle->GetKineticEnergy();
82  G4cout<<"G4QIonIonElastic::GetMeanFreePath: kinE="<<KinEn<<",Mom="<<Momentum<<G4endl;
83#endif
84  const G4Material* material = aTrack.GetMaterial();        // Get the current material
85  const G4double* NOfNucPerVolume = material->GetVecNbOfAtomsPerVolume();
86  const G4ElementVector* theElementVector = material->GetElementVector();
87  G4int nE=material->GetNumberOfElements();
88#ifdef debug
89  G4cout<<"G4QIonIonElastic::GetMeanFreePath:"<<nE<<" Elem's in theMaterial"<<G4endl;
90#endif
91  G4VQCrossSection* CSmanager=G4QIonIonCrossSection::GetPointer();
92  G4int pPDG=0;
93  // Probably enough: pPDG=incidentParticleDefinition->GetPDGEncoding();
94  if      ( incidentParticleDefinition ==  G4Deuteron::Deuteron()     ) pPDG = 100001002;
95  else if ( incidentParticleDefinition ==  G4Alpha::Alpha()           ) pPDG = 100002004;
96  else if ( incidentParticleDefinition ==  G4Triton::Triton()         ) pPDG = 100001003;
97  else if ( incidentParticleDefinition ==  G4He3::He3()               ) pPDG = 100002003;
98  else if ( incidentParticleDefinition ==  G4GenericIon::GenericIon() )
99  {
100    pPDG=incidentParticleDefinition->GetPDGEncoding();
101#ifdef debug
102    G4int B=incidentParticleDefinition->GetBaryonNumber();
103    G4int C=incidentParticleDefinition->GetPDGCharge();
104    prPDG=100000000+1000*C+B;
105    G4cout<<"G4QIonIonElastic::GetMeanFreePath: PDG="<<prPDG<<"="<<pPDG<<G4endl;
106#endif
107  }
108  else G4cout<<"-Warning-G4QIonIonElastic::GetMeanFreePath:Unknown projectile Ion"<<G4endl;
109  Momentum/=incidentParticleDefinition->GetBaryonNumber(); // Divide Mom by projectile A
110  G4QIsotope* Isotopes = G4QIsotope::Get(); // Pointer to the G4QIsotopes singleton
111  G4double sigma=0.;                        // Sums over elements for the material
112  G4int IPIE=IsoProbInEl.size();            // How many old elements?
113  if(IPIE) for(G4int ip=0; ip<IPIE; ++ip)   // Clean up the SumProb's of Isotopes (SPI)
114  {
115    std::vector<G4double>* SPI=IsoProbInEl[ip]; // Pointer to the SPI vector
116    SPI->clear();
117    delete SPI;
118    std::vector<G4int>* IsN=ElIsoN[ip];     // Pointer to the N vector
119    IsN->clear();
120    delete IsN;
121  }
122  ElProbInMat.clear();                      // Clean up the SumProb's of Elements (SPE)
123  ElementZ.clear();                         // Clear the body vector for Z of Elements
124  IsoProbInEl.clear();                      // Clear the body vector for SPI
125  ElIsoN.clear();                           // Clear the body vector for N of Isotopes
126  for(G4int i=0; i<nE; ++i)
127  {
128    G4Element* pElement=(*theElementVector)[i]; // Pointer to the current element
129    G4int Z = static_cast<G4int>(pElement->GetZ()); // Z of the Element
130    ElementZ.push_back(Z);                  // Remember Z of the Element
131    G4int isoSize=0;                        // The default for the isoVectorLength is 0
132    G4int indEl=0;                          // Index of non-natural element or 0(default)
133    G4IsotopeVector* isoVector=pElement->GetIsotopeVector(); // Get the predefined IsoVect
134    if(isoVector) isoSize=isoVector->size();// Get size of the existing isotopeVector
135#ifdef debug
136    G4cout<<"G4QIonIonElastic::GetMeanFreePath: isovectorLength="<<isoSize<<G4endl;
137#endif
138    if(isoSize)                             // The Element has non-trivial abundance set
139    {
140      indEl=pElement->GetIndex()+1;         // Index of the non-trivial element is an order
141#ifdef debug
142      G4cout<<"G4QIIEl::GetMFP:iE="<<indEl<<", def="<<Isotopes->IsDefined(Z,indEl)<<G4endl;
143#endif
144      if(!Isotopes->IsDefined(Z,indEl))     // This index is not defined for this Z: define
145      {
146        std::vector<std::pair<G4int,G4double>*>* newAbund =
147                                               new std::vector<std::pair<G4int,G4double>*>;
148        G4double* abuVector=pElement->GetRelativeAbundanceVector();
149        for(G4int j=0; j<isoSize; j++)      // Calculation of abundance vector for isotopes
150        {
151          G4int N=pElement->GetIsotope(j)->GetN()-Z; // N means A=N+Z !
152          if(pElement->GetIsotope(j)->GetZ()!=Z) G4cerr<<"G4QIonIonEl::GetMeanFreePath Z="
153                                         <<pElement->GetIsotope(j)->GetZ()<<"#"<<Z<<G4endl;
154          G4double abund=abuVector[j];
155                                                                  std::pair<G4int,G4double>* pr= new std::pair<G4int,G4double>(N,abund);
156#ifdef debug
157          G4cout<<"G4QIonIonElastic::GetMeanFP:pair#="<<j<<",N="<<N<<",ab="<<abund<<G4endl;
158#endif
159          newAbund->push_back(pr);
160                                                  }
161#ifdef debug
162        G4cout<<"G4QIonIonElastic::GetMeanFP: pairVectorLength="<<newAbund->size()<<G4endl;
163#endif
164        indEl=G4QIsotope::Get()->InitElement(Z,indEl,newAbund); // definition of the newInd
165        for(G4int k=0; k<isoSize; k++) delete (*newAbund)[k];   // Cleaning temporary
166        delete newAbund; // Was "new" in the beginning of the name space
167      }
168    }
169    std::vector<std::pair<G4int,G4double>*>* cs= Isotopes->GetCSVector(Z,indEl);//CSPointer
170    std::vector<G4double>* SPI = new std::vector<G4double>; // Pointer to the SPI vector
171    IsoProbInEl.push_back(SPI);
172    std::vector<G4int>* IsN = new std::vector<G4int>; // Pointer to the N vector
173    ElIsoN.push_back(IsN);
174    G4int nIs=cs->size();                   // A#Of Isotopes in the Element
175#ifdef debug
176    G4cout<<"G4QIonIonEl::GetMFP:=***=> #isot="<<nIs<<", Z="<<Z<<", indEl="<<indEl<<G4endl;
177#endif
178    G4double susi=0.;                       // sum of CS over isotopes
179    if(nIs) for(G4int j=0; j<nIs; j++)      // Calculate CS for eachIsotope of El
180    {
181      std::pair<G4int,G4double>* curIs=(*cs)[j]; // A pointer, which is used twice
182      G4int N=curIs->first;                 // #of Neuterons in the isotope j of El i
183      IsN->push_back(N);                    // Remember Min N for the Element
184#ifdef debug
185      G4cout<<"G4QIIEl::GMFP:true,P="<<Momentum<<",Z="<<Z<<",N="<<N<<",PDG="<<pPDG<<G4endl;
186#endif
187                    G4bool ccsf=false;                    // Extract elastic Ion-Ion cross-section
188#ifdef debug
189      G4cout<<"G4QIonIonElastic::GMFP: GetCS #1 j="<<j<<G4endl;
190#endif
191      G4double CSI=CSmanager->GetCrossSection(ccsf,Momentum,Z,N,pPDG);//CS(j,i) for isotope
192#ifdef debug
193      G4cout<<"G4QIonIonElastic::GMFP: jI="<<j<<", Zt="<<Z<<", Nt="<<N<<", Mom="<<Momentum
194            <<", XSec="<<CSI/millibarn<<G4endl;
195#endif
196      curIs->second = CSI;
197      susi+=CSI;                            // Make a sum per isotopes
198      SPI->push_back(susi);                 // Remember summed cross-section
199    } // End of temporary initialization of the cross sections in the G4QIsotope singeltone
200    sigma+=Isotopes->GetMeanCrossSection(Z,indEl)*NOfNucPerVolume[i];//SUM(MeanCS*NOfNperV)
201#ifdef debug
202    G4cout<<"G4QIonIonEl::GMFP:<S>="<<Isotopes->GetMeanCrossSection(Z,indEl)<<",AddToSig="
203          <<Isotopes->GetMeanCrossSection(Z,indEl)*NOfNucPerVolume[i]<<G4endl;
204#endif
205    ElProbInMat.push_back(sigma);
206  } // End of LOOP over Elements
207  // Check that cross section is not zero and return the mean free path
208#ifdef debug
209  G4cout<<"G4QIonIonElastic::GetMeanFreePath: MeanFreePath="<<1./sigma<<G4endl;
210#endif
211  if(sigma > 0.) return 1./sigma;                 // Mean path [distance]
212  return DBL_MAX;
213}
214
215
216G4bool G4QIonIonElastic::IsApplicable(const G4ParticleDefinition& particle) 
217{
218  if      (particle == *( G4Deuteron::Deuteron()     )) return true;
219  else if (particle == *( G4Alpha::Alpha()           )) return true;
220  else if (particle == *( G4Triton::Triton()         )) return true;
221  else if (particle == *( G4He3::He3()               )) return true;
222  else if (particle == *( G4GenericIon::GenericIon() )) return true;
223#ifdef debug
224  G4cout<<"***>>G4QIonIonElastic::IsApplicable: PDG="<<particle.GetPDGEncoding()<<G4endl;
225#endif
226  return false;
227}
228
229G4VParticleChange* G4QIonIonElastic::PostStepDoIt(const G4Track& track, const G4Step& step)
230{
231  static const G4double fm2MeV2 = 3*38938./1.09; // (3/1.09)*(hc)^2 in fm^2*MeV^2
232  static G4bool CWinit = true;                   // CHIPS Warld needs to be initted
233  if(CWinit)
234                {
235    CWinit=false;
236    G4QCHIPSWorld::Get()->GetParticles(nPartCWorld); // Create CHIPS World (234 part.max)
237  }
238  //-------------------------------------------------------------------------------------
239  const G4DynamicParticle* projHadron = track.GetDynamicParticle();
240  const G4ParticleDefinition* particle=projHadron->GetDefinition();
241#ifdef debug
242  G4cout<<"G4QIonIonElastic::PostStepDoIt: Before the GetMeanFreePath is called In4M="
243        <<projHadron->Get4Momentum()<<" of PDG="<<particle->GetPDGEncoding()<<", Type="
244        <<particle->GetParticleType()<<", Subtp="<<particle->GetParticleSubType()<<G4endl;
245#endif
246  //G4ForceCondition cond=NotForced;
247  //GetMeanFreePath(track, -27., &cond);                // @@ ?? jus to update parameters?
248#ifdef debug
249  G4cout<<"G4QIonIonElastic::PostStepDoIt: After the GetMeanFreePath is called"<<G4endl;
250#endif
251  G4LorentzVector proj4M=(projHadron->Get4Momentum())/MeV; // Convert to MeV!
252  G4LorentzVector scat4M=proj4M;                      // @@ Must be filled (?)
253  G4double momentum = projHadron->GetTotalMomentum()/MeV; // 3-momentum of the Proj in MeV
254  G4double Momentum = proj4M.rho();                   // @@ Just for the test purposes
255  if(std::fabs(Momentum-momentum)>.000001)
256    G4cerr<<"-Warn-G4QIonIonElastic::PostStepDoIt:P(IU)="<<Momentum<<"#"<<momentum<<G4endl;
257  G4double pM2=proj4M.m2();        // in MeV^2
258  G4double pM=std::sqrt(pM2);      // in MeV
259#ifdef pdebug
260  G4cout<<"G4QIonIonElastic::PostStepDoIt: pP(IU)="<<Momentum<<"="<<momentum<<",pM="<<pM
261        <<",scat4M="<<scat4M<<scat4M.m()<<G4endl;
262#endif
263  if (!IsApplicable(*particle))  // Check applicability
264  {
265    G4cerr<<"G4QIonIonElastic::PostStepDoIt: Only NA elastic is implemented."<<G4endl;
266    return 0;
267  }
268  const G4Material* material = track.GetMaterial();      // Get the current material
269  G4int Z=0;
270  const G4ElementVector* theElementVector = material->GetElementVector();
271  G4int nE=material->GetNumberOfElements();
272#ifdef debug
273  G4cout<<"G4QIonIonElastic::PostStepDoIt: "<<nE<<" elements in the material."<<G4endl;
274#endif
275  // Probably enough: projPDG=particle->GetPDGEncoding();
276  G4int projPDG=0;                           // CHIPS PDG Code for the captured hadron
277  if      (particle ==  G4Deuteron::Deuteron()     ) projPDG= 100001002;
278  else if (particle ==  G4Alpha::Alpha()           ) projPDG= 100002004;
279  else if (particle ==  G4Triton::Triton()         ) projPDG= 100001003;
280  else if (particle ==  G4He3::He3()               ) projPDG= 100002003;
281  else if (particle ==  G4GenericIon::GenericIon() )
282  {
283    projPDG=particle->GetPDGEncoding();
284#ifdef debug
285    G4int B=particle->GetBaryonNumber();
286    G4int C=particle->GetPDGCharge();
287    prPDG=100000000+1000*C+B;
288    G4cout<<"G4QIonIonElastic::PostStepDoIt: PDG="<<prPDG<<"="<<projPDG<<G4endl;
289#endif
290  }
291  else G4cout<<"-Warning-G4QIonIonElastic::PostStepDoIt:Unknown projectile Ion"<<G4endl;
292#ifdef debug
293  G4int prPDG=particle->GetPDGEncoding();
294                G4cout<<"G4QIonIonElastic::PostStepDoIt: projPDG="<<projPDG<<", stPDG="<<prPDG<<G4endl;
295#endif
296  if(!projPDG)
297  {
298    G4cerr<<"-Warning-G4QIonIonElastic::PostStepDoIt:Undefined interactingNucleus"<<G4endl;
299    return 0;
300  }
301  G4double pA=particle->GetBaryonNumber();     // Projectile A
302  G4double pZ=particle->GetPDGCharge();        // Projectile Z
303  G4double pN=pA-pZ;                           // Projectile N
304  G4int EPIM=ElProbInMat.size();
305#ifdef debug
306                G4cout<<"G4QIonIonElastic::PSDI:m="<<EPIM<<",n="<<nE<<",T="<<ElProbInMat[EPIM-1]<<G4endl;
307#endif
308  G4int i=0;
309  if(EPIM>1)
310  {
311    G4double rnd = ElProbInMat[EPIM-1]*G4UniformRand();
312    for(i=0; i<nE; ++i)
313                  {
314#ifdef debug
315                                  G4cout<<"G4QIonIonElastic::PSDI: EPM["<<i<<"]="<<ElProbInMat[i]<<", r="<<rnd<<G4endl;
316#endif
317      if (rnd<ElProbInMat[i]) break;
318    }
319    if(i>=nE) i=nE-1;                        // Top limit for the Element
320  }
321  G4Element* pElement=(*theElementVector)[i];
322  Z=static_cast<G4int>(pElement->GetZ());
323#ifdef debug
324                                G4cout<<"G4QIonIonElastic::PostStepDoIt: i="<<i<<", Z(element)="<<Z<<G4endl;
325#endif
326  if(Z<=0)
327  {
328    G4cerr<<"---Warning---G4QIonIonElastic::PostStepDoIt: Element with Z="<<Z<<G4endl;
329    if(Z<0) return 0;
330  }
331  std::vector<G4double>* SPI = IsoProbInEl[i];// Vector of summedProbabilities for isotopes
332  std::vector<G4int>* IsN = ElIsoN[i];     // Vector of "#of neutrons" in the isotope El[i]
333  G4int nofIsot=SPI->size();               // #of isotopes in the element i
334#ifdef debug
335                G4cout<<"G4QIonIonElastic::PosStDoIt: nI="<<nofIsot<<",T="<<(*SPI)[nofIsot-1]<<G4endl;
336#endif
337  G4int j=0;
338  if(nofIsot>1)
339  {
340    G4double rndI=(*SPI)[nofIsot-1]*G4UniformRand(); // Randomize the isotop of the Element
341    for(j=0; j<nofIsot; ++j)
342    {
343#ifdef debug
344                                  G4cout<<"G4QIonIonElastic::PostStDI: SP["<<j<<"]="<<(*SPI)[j]<<", r="<<rndI<<G4endl;
345#endif
346      if(rndI < (*SPI)[j]) break;
347    }
348    if(j>=nofIsot) j=nofIsot-1;            // Top limit for the isotope
349  }
350  G4int N =(*IsN)[j]; ;                    // Randomized number of neutrons
351#ifdef debug
352                G4cout<<"G4QIonIonElastic::PostStepDoIt:j="<<i<<",N(isotope)="<<N<<", MeV="<<MeV<<G4endl;
353#endif
354  if(N<0)
355  {
356    G4cerr<<"-Warning-G4QIonIonElastic::PostStepDoIt:IsotopeZ="<<Z<<" has 0>N="<<N<<G4endl;
357    return 0;
358  }
359  nOfNeutrons=N;                           // Remember it for the energy-momentum check
360#ifdef debug
361  G4cout<<"G4QIonIonElastic::PostStepDoIt: N="<<N<<" for element with Z="<<Z<<G4endl;
362#endif
363  if(N<0)
364  {
365    G4cerr<<"-Warning-G4QIonIonElastic::PostStepDoIt:Element with N="<<N<< G4endl;
366    return 0;
367  }
368  aParticleChange.Initialize(track);
369#ifdef debug
370  G4cout<<"G4QIonIonElastic::PostStepDoIt: track is initialized"<<G4endl;
371#endif
372  G4double weight        = track.GetWeight();
373  G4double localtime     = track.GetGlobalTime();
374  G4ThreeVector position = track.GetPosition();
375#ifdef debug
376  G4cout<<"G4QIonIonElastic::PostStepDoIt: before Touchable extraction"<<G4endl;
377#endif
378  G4TouchableHandle trTouchable = track.GetTouchableHandle();
379#ifdef debug
380  G4cout<<"G4QIonIonElastic::PostStepDoIt: Touchable is extracted"<<G4endl;
381#endif
382  //
383  G4double tA=Z+N;
384  G4int targPDG=90000000+Z*1000+N;         // CHIPS PDG Code of the target nucleus
385  G4QPDGCode targQPDG(targPDG);            // @@ one can use G4Ion & get rid of CHIPS World
386  G4double tM=targQPDG.GetMass();          // CHIPS target mass in MeV
387  G4double kinEnergy= projHadron->GetKineticEnergy()*MeV; // Kin energy in MeV (Is *MeV n?)
388  G4ParticleMomentum dir = projHadron->GetMomentumDirection();// It is a unit three-vector
389  G4LorentzVector tot4M=proj4M+G4LorentzVector(0.,0.,0.,tM); // Total 4-mom of the reaction
390#ifdef debug
391  G4cout<<"G4QIonIonElastic::PostStDI: tM="<<tM<<", p4M="<<proj4M<<", t4M="<<tot4M<<G4endl;
392#endif
393  EnMomConservation=tot4M;                 // Total 4-mom of reaction for E/M conservation
394  G4VQCrossSection* ELmanager=G4QElasticCrossSection::GetPointer();
395  G4VQCrossSection* CSmanager=G4QIonIonCrossSection::GetPointer();
396  // @@ Probably this is not necessary any more
397#ifdef debug
398  G4cout<<"G4QIIEl::PSDI:false,P="<<Momentum<<",Z="<<Z<<",N="<<N<<",PDG="<<projPDG<<G4endl;
399#endif
400  // false means elastic cross-section
401  G4double xSec=CSmanager->GetCrossSection(false, Momentum, Z, N, projPDG);// Rec.CrossSect
402#ifdef debug
403  G4cout<<"G4QIIEl::PSDI: pPDG="<<projPDG<<",P="<<Momentum<<",CS="<<xSec/millibarn<<G4endl;
404#endif
405#ifdef nandebug
406  if(xSec>0. || xSec<0. || xSec==0);
407  else  G4cout<<"-NaN-Warning-G4QIonIonElastic::PostStDoIt: xSec="<<xSec/millibarn<<G4endl;
408#endif
409  // @@ check a possibility to separate p, n, or alpha (!)
410  if(xSec <= 0.) // The cross-section iz 0 -> Do Nothing
411  {
412#ifdef pdebug
413    G4cerr<<"-Warning-G4QIonIonElastic::PostStDoIt: *Zero cross-section* PDG="<<projPDG
414          <<",tPDG="<<targPDG<<",P="<<Momentum<<G4endl;
415#endif
416    //Do Nothing Action insead of the reaction
417    aParticleChange.ProposeEnergy(kinEnergy);
418    aParticleChange.ProposeLocalEnergyDeposit(0.);
419    aParticleChange.ProposeMomentumDirection(dir) ;
420    return G4VDiscreteProcess::PostStepDoIt(track,step);
421  }
422  G4double dtM=tM+tM;
423  G4double PA=Momentum*pA;
424  G4double PA2=PA*PA;
425  G4double maxt=dtM*PA2/(std::sqrt(PA2+pM2)+tM/2+pM2/dtM);
426#ifdef pdebug
427  G4cout<<"G4QIonIonElastic::PostStDoIt:pPDG="<<projPDG<<",tPDG="<<targPDG<<",P="<<Momentum
428        <<",CS="<<xSec<<",maxt="<<maxt<<G4endl;
429#endif
430  xSec=ELmanager->GetCrossSection(false, Momentum, 1, 0, 2212);// pp=nn
431  G4double B1=ELmanager->GetSlope(1,0,2212); // slope for pp=nn
432  xSec=ELmanager->GetCrossSection(false, Momentum, 1, 0, 2112);// np=pn
433  G4double B2 =ELmanager->GetSlope(1,0,2112); // slope for np=pn
434  G4double mB =((pZ*Z+pN*N)*B1+(pZ*N+pN*Z)*B2)/(pA+tA);
435  G4double pR2=std::pow(pA+4.,.305)/fm2MeV2;
436  G4double tR2=std::pow(tA+4.,.305)/fm2MeV2;
437  G4double eB =mB+pR2+tR2;
438  G4double mint=-std::log(1.-G4UniformRand()*(1.-std::exp(-eB*maxt)))/eB;
439#ifdef pdebug
440  G4cout<<"G4QIonIonElastic::PostStDoIt:B1="<<B1<<",B2="<<B2<<",mB="<<mB
441        <<",pR2="<<pR2<<",tR2="<<tR2<<",eB="<<eB<<",mint="<<mint<<G4endl;
442#endif
443#ifdef nandebug
444  if(mint>-.0000001);
445  else  G4cout<<"-Warning-G4QIonIonElastic::PostStDoIt:-t="<<mint<<G4endl;
446#endif
447  G4double cost=1.-(mint+mint)/maxt; // cos(theta) in CMS
448  //
449#ifdef ppdebug
450  G4cout<<"G4QIonIonElastic::PoStDoI:t="<<mint<<",dpcm2="<<CSmanager->GetHMaxT()<<",Ek="
451        <<kinEnergy<<",tM="<<tM<<",pM="<<pM<<",cost="<<cost<<G4endl;
452#endif
453  if(cost>1. || cost<-1. || !(cost>-1. || cost<1.))
454  {
455    if(cost>1.000001 || cost<-1.000001 || !(cost>-1. || cost<1.))
456    {
457      G4double tM2=tM*tM;                         // Squared target mass
458      G4double pEn=pM+kinEnergy;                  // tot projectile Energy in MeV
459      G4double sM=dtM*pEn+tM2+pM2;                // Mondelstam s
460      G4double twop2cm=(tM2+tM2)*(pEn*pEn-pM2)/sM;// Max_t/2 (2*p^2_cm)
461      G4cout<<"-Warning-G4QIonIonElastic::PoStDI:cos="<<cost<<",t="<<mint<<",T="<<kinEnergy
462            <<",tM="<<tM<<",tmax="<<2*kinEnergy*tM<<",p="<<projPDG<<",t="<<targPDG<<G4endl;
463      G4cout<<"G4QIonIonElastic::PSDI:dpcm2="<<twop2cm<<"="<<CSmanager->GetHMaxT()<<G4endl;
464    }
465    if     (cost>1.)  cost=1.;
466    else if(cost<-1.) cost=-1.;
467  }
468  G4LorentzVector reco4M=G4LorentzVector(0.,0.,0.,tM);      // 4mom of the recoil target
469  G4LorentzVector dir4M=tot4M-G4LorentzVector(0.,0.,0.,(tot4M.e()-tM-pM)*.01);
470  if(!G4QHadron(tot4M).RelDecayIn2(scat4M, reco4M, dir4M, cost, cost))
471  {
472    G4cerr<<"G4QIonIonE::PSDI:t4M="<<tot4M<<",pM="<<pM<<",tM="<<tM<<",cost="<<cost<<G4endl;
473  }
474#ifdef debug
475  G4cout<<"G4QIonIonElast::PSDI:s4M="<<scat4M<<"+r4M="<<reco4M<<"="<<scat4M+reco4M<<G4endl;
476  G4cout<<"G4QIonIonElastic::PSDI:scatE="<<scat4M.e()-pM<<",recoE="<<reco4M.e()-tM<<",d4M="
477        <<tot4M-scat4M-reco4M<<G4endl;
478#endif
479  // Update G4VParticleChange for the scattered projectile
480  G4double finE=scat4M.e()-pM;             // Final kinetic energy of the scattered proton
481  if(finE>0.0) aParticleChange.ProposeEnergy(finE);
482  else
483  {
484    if(finE<-1.e-8 || !(finE>-1.||finE<1.)) // NAN or negative
485      G4cerr<<"*Warning*G4QIonIonElastic::PostStDoIt: Zero or negative scattered E="<<finE
486            <<", s4M="<<scat4M<<", r4M="<<reco4M<<", d4M="<<tot4M-scat4M-reco4M<<G4endl;
487    aParticleChange.ProposeEnergy(0.) ;
488    aParticleChange.ProposeTrackStatus(fStopButAlive);
489  }
490  G4ThreeVector findir=scat4M.vect()/scat4M.rho();  // Unit vector in new direction
491  aParticleChange.ProposeMomentumDirection(findir); // new direction for the scattered part
492  EnMomConservation-=scat4M;                        // It must be initialized by (pE+tM,pP)
493  // This is how in general the secondary should be identified
494                G4DynamicParticle* theSec = new G4DynamicParticle; // A secondary for the recoil hadron
495  G4int aA = Z+N;
496#ifdef pdebug
497                G4cout<<"G4QIonIonElastic::PostStepDoIt: Ion Z="<<Z<<", A="<<aA<<G4endl;
498#endif
499  G4ParticleDefinition* theDefinition=G4ParticleTable::GetParticleTable()
500                                                                       ->FindIon(Z,aA,0,Z);
501  if(!theDefinition)G4cout<<"-Warning-G4QIonIonElastic::PoStDI:drop PDG="<<targPDG<<G4endl;
502#ifdef pdebug
503  G4cout<<"G4QIonIonElastic::PoStDI:RecoilName="<<theDefinition->GetParticleName()<<G4endl;
504#endif
505  theSec->SetDefinition(theDefinition);
506
507  EnMomConservation-=reco4M;
508#ifdef tdebug
509  G4cout<<"G4QIonIonElastic::PSD:"<<targPDG<<reco4M<<reco4M.m()<<EnMomConservation<<G4endl;
510#endif
511  theSec->Set4Momentum(reco4M);
512#ifdef debug
513  G4ThreeVector curD=theSec->GetMomentumDirection();
514  G4double curM=theSec->GetMass();
515  G4double curE=theSec->GetKineticEnergy()+curM;
516  G4cout<<"G4QIonIonElastic::PSDI: p="<<curD<<curD.mag()<<",e="<<curE<<",m="<<curM<<G4endl;
517#endif
518  // Make a recoil nucleus
519  G4Track* aNewTrack = new G4Track(theSec, localtime, position );
520  aNewTrack->SetWeight(weight);                                   //    weighted
521  aNewTrack->SetTouchableHandle(trTouchable);
522  aParticleChange.AddSecondary( aNewTrack );
523#ifdef debug
524    G4cout<<"G4QIonIonElastic::PostStepDoIt: **** PostStepDoIt is done ****"<<G4endl;
525#endif
526  return G4VDiscreteProcess::PostStepDoIt(track, step);
527}
Note: See TracBrowser for help on using the repository browser.