source: trunk/source/processes/hadronic/models/chiral_inv_phase_space/interface/src/G4QCoherentChargeExchange.cc @ 1197

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

maj sur la beta de geant 4.9.3

File size: 30.0 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: G4QCoherentChargeExchange.cc,v 1.8 2009/02/23 09:49:24 mkossov Exp $
27// GEANT4 tag $Name: geant4-09-03-beta-cand-01 $
28//
29//      ---------------- G4QCoherentChargeExchange class -----------------
30//                 by Mikhail Kossov, December 2003.
31// G4QCoherentChargeExchange class of the CHIPS Simulation Branch in GEANT4
32// ---------------------------------------------------------------
33// ****************************************************************************************
34// ********** This CLASS is temporary moved from the photolepton_hadron directory *********
35// ****************************************************************************************
36// Short description: This class resolves an ambiguity in the definition of the
37// "inelastic" cross section. As it was shown in Ph.D.Thesis (M.Kosov,ITEP,1979)
38// it is more reasonable to subdivide the total cross-section in the coherent &
39// incoherent parts, but the measuring method for the "inelastic" cross-sections
40// consideres the lack of the projectile within the narrow forward solid angle
41// with the consequent extrapolation of these partial cross-sections, corresponding
42// to the particular solid angle, to the zero solid angle. The low angle region
43// is shadowed by the elastic (coherent) scattering. BUT the coherent charge
44// exchange (e.g. conversion p->n) is included by this procedure as a constant term
45// in the extrapolation, so the "inelastic" cross-section differes from the
46// incoherent cross-section by the value of the coherent charge exchange cross
47// section. Fortunately, this cross-sectoion drops ruther fast with energy increasing.
48// All Geant4 inelastic hadronic models (including CHIPS) simulate the incoherent
49// reactions. So the incoherent (including quasielastic) cross-section must be used
50// instead of the inelastic cross-section. For that the "inelastic" cross-section
51// must be reduced by the value of the coherent charge-exchange cross-section, which
52// is estimated (it must be tuned!) in this CHIPS class. The angular distribution
53// is made (at present) identical to the corresponding coherent-elastic scattering
54// -----------------------------------------------------------------------------------
55//#define debug
56//#define pdebug
57//#define tdebug
58//#define nandebug
59//#define ppdebug
60
61#include "G4QCoherentChargeExchange.hh"
62
63// Initialization of static vectors
64G4int    G4QCoherentChargeExchange::nPartCWorld=152;  // #of particles initialized in CHIPS
65std::vector<G4int> G4QCoherentChargeExchange::ElementZ; // Z of element(i) in theLastCalc
66std::vector<G4double> G4QCoherentChargeExchange::ElProbInMat; // SumProbOfElem in Material
67std::vector<std::vector<G4int>*> G4QCoherentChargeExchange::ElIsoN;// N of isotope(j), E(i)
68std::vector<std::vector<G4double>*>G4QCoherentChargeExchange::IsoProbInEl;//SumProbIsotE(i)
69
70// Constructor
71G4QCoherentChargeExchange::G4QCoherentChargeExchange(const G4String& processName)
72  : G4VDiscreteProcess(processName, fHadronic)
73{
74#ifdef debug
75  G4cout<<"G4QCohChargeEx::Constructor is called processName="<<processName<<G4endl;
76#endif
77  if (verboseLevel>0) G4cout << GetProcessName() << " process is created "<< G4endl;
78  SetProcessSubType(fChargeExchange);
79  //G4QCHIPSWorld::Get()->GetParticles(nPartCWorld); // Create CHIPS World (234 part. max)
80}
81
82// Destructor
83G4QCoherentChargeExchange::~G4QCoherentChargeExchange() {}
84
85
86G4LorentzVector G4QCoherentChargeExchange::GetEnegryMomentumConservation()
87                                                                {return EnMomConservation;}
88
89G4int G4QCoherentChargeExchange::GetNumberOfNeutronsInTarget() {return nOfNeutrons;}
90
91// output of the function must be in units of length! L=1/sig_V,sig_V=SUM(n(j,i)*sig(j,i)),
92// where n(i,j) is a number of nuclei of the isotop j of the element i in V=1(lengtUnit^3)
93// ********** All CHIPS cross sections are calculated in the surface units ************
94G4double G4QCoherentChargeExchange::GetMeanFreePath(const G4Track& aTrack,G4double Q,
95                                                    G4ForceCondition* Fc)
96{
97  *Fc = NotForced;
98  const G4DynamicParticle* incidentParticle = aTrack.GetDynamicParticle();
99  G4ParticleDefinition* incidentParticleDefinition=incidentParticle->GetDefinition();
100  if( !IsApplicable(*incidentParticleDefinition))
101    G4cout<<"*W*G4QCohChargeEx::GetMeanFreePath called for notImplementedParticle"<<G4endl;
102  // Calculate the mean Cross Section for the set of Elements(*Isotopes) in the Material
103  G4double Momentum = incidentParticle->GetTotalMomentum(); // 3-momentum of the Particle
104#ifdef debug
105  G4double KinEn = incidentParticle->GetKineticEnergy();
106  G4cout<<"G4QCohChEx::GetMeanFreePath: kinE="<<KinEn<<",Mom="<<Momentum<<G4endl; // Result
107#endif
108  const G4Material* material = aTrack.GetMaterial();        // Get the current material
109  const G4double* NOfNucPerVolume = material->GetVecNbOfAtomsPerVolume();
110  const G4ElementVector* theElementVector = material->GetElementVector();
111  G4int nE=material->GetNumberOfElements();
112#ifdef debug
113  G4cout<<"G4QCohChargeExchange::GetMeanFreePath:"<<nE<<" Elem's in theMaterial"<<G4endl;
114#endif
115  G4int pPDG=0;
116
117  if     (incidentParticleDefinition == G4Proton::Proton()  ) pPDG=2212;
118  else if(incidentParticleDefinition == G4Neutron::Neutron()) pPDG=2112;
119  else G4cout<<"G4QCohChargeEx::GetMeanFreePath: only nA & pA are implemented"<<G4endl;
120 
121  G4QIsotope* Isotopes = G4QIsotope::Get(); // Pointer to the G4QIsotopes singleton
122  G4double sigma=0.;                        // Sums over elements for the material
123  G4int IPIE=IsoProbInEl.size();            // How many old elements?
124  if(IPIE) for(G4int ip=0; ip<IPIE; ++ip)   // Clean up the SumProb's of Isotopes (SPI)
125  {
126    std::vector<G4double>* SPI=IsoProbInEl[ip]; // Pointer to the SPI vector
127    SPI->clear();
128    delete SPI;
129    std::vector<G4int>* IsN=ElIsoN[ip];     // Pointer to the N vector
130    IsN->clear();
131    delete IsN;
132  }
133  ElProbInMat.clear();                      // Clean up the SumProb's of Elements (SPE)
134  ElementZ.clear();                         // Clear the body vector for Z of Elements
135  IsoProbInEl.clear();                      // Clear the body vector for SPI
136  ElIsoN.clear();                           // Clear the body vector for N of Isotopes
137  for(G4int i=0; i<nE; ++i)
138  {
139    G4Element* pElement=(*theElementVector)[i]; // Pointer to the current element
140    G4int Z = static_cast<G4int>(pElement->GetZ()); // Z of the Element
141    ElementZ.push_back(Z);                  // Remember Z of the Element
142    G4int isoSize=0;                        // The default for the isoVectorLength is 0
143    G4int indEl=0;                          // Index of non-natural element or 0(default)
144    G4IsotopeVector* isoVector=pElement->GetIsotopeVector(); // Get the predefined IsoVect
145    if(isoVector) isoSize=isoVector->size();// Get size of the existing isotopeVector
146#ifdef debug
147    G4cout<<"G4QCoherentChargeExchange::GetMeanFreePath:isovectorLength="<<isoSize<<G4endl;
148#endif
149    if(isoSize)                             // The Element has non-trivial abundance set
150    {
151      indEl=pElement->GetIndex()+1;         // Index of the non-trivial element is an order
152#ifdef debug
153      G4cout<<"G4QCCX::GetMFP: iE="<<indEl<<", def="<<Isotopes->IsDefined(Z,indEl)<<G4endl;
154#endif
155      if(!Isotopes->IsDefined(Z,indEl))     // This index is not defined for this Z: define
156      {
157        std::vector<std::pair<G4int,G4double>*>* newAbund =
158                                               new std::vector<std::pair<G4int,G4double>*>;
159        G4double* abuVector=pElement->GetRelativeAbundanceVector();
160        for(G4int j=0; j<isoSize; j++)      // Calculation of abundance vector for isotopes
161        {
162          G4int N=pElement->GetIsotope(j)->GetN()-Z; // N means A=N+Z !
163          if(pElement->GetIsotope(j)->GetZ()!=Z)G4cerr<<"G4QCohChX::GetMeanFreePath: Z="
164                                         <<pElement->GetIsotope(j)->GetZ()<<"#"<<Z<<G4endl;
165          G4double abund=abuVector[j];
166          std::pair<G4int,G4double>* pr= new std::pair<G4int,G4double>(N,abund);
167#ifdef debug
168          G4cout<<"G4QCohChEx::GetMeanFreePath:pair#="<<j<<",N="<<N<<",ab="<<abund<<G4endl;
169#endif
170          newAbund->push_back(pr);
171        }
172#ifdef debug
173        G4cout<<"G4QCohChEx::GetMeanFreePath: pairVectorLength="<<newAbund->size()<<G4endl;
174#endif
175        indEl=G4QIsotope::Get()->InitElement(Z,indEl,newAbund); // definition of the newInd
176        for(G4int k=0; k<isoSize; k++) delete (*newAbund)[k];   // Cleaning temporary
177        delete newAbund; // Was "new" in the beginning of the name space
178      }
179    }
180    std::vector<std::pair<G4int,G4double>*>* cs= Isotopes->GetCSVector(Z,indEl);//CSPointer
181    std::vector<G4double>* SPI = new std::vector<G4double>; // Pointer to the SPI vector
182    IsoProbInEl.push_back(SPI);
183    std::vector<G4int>* IsN = new std::vector<G4int>; // Pointer to the N vector
184    ElIsoN.push_back(IsN);
185    G4int nIs=cs->size();                   // A#Of Isotopes in the Element
186#ifdef debug
187    G4cout<<"G4QCohChargEx::GMFP:=***=>,#isot="<<nIs<<", Z="<<Z<<", indEl="<<indEl<<G4endl;
188#endif
189    G4double susi=0.;                       // sum of CS over isotopes
190    if(nIs) for(G4int j=0; j<nIs; j++)      // Calculate CS for eachIsotope of El
191    {
192      std::pair<G4int,G4double>* curIs=(*cs)[j]; // A pointer, which is used twice
193      G4int N=curIs->first;                 // #of Neuterons in the isotope j of El i
194      IsN->push_back(N);                    // Remember Min N for the Element
195#ifdef debug
196      G4cout<<"G4QCCX::GMFP:true, P="<<Momentum<<",Z="<<Z<<",N="<<N<<",PDG="<<pPDG<<G4endl;
197#endif
198      G4bool ccsf=true;
199      if(Q==-27.) ccsf=false;
200#ifdef debug
201      G4cout<<"G4QCoherentChargeExchange::GMFP: GetCS #1 j="<<j<<G4endl;
202#endif
203      G4double CSI=CalculateXSt(ccsf, true, Momentum, Z, N, pPDG);// CS(j,i) for theIsotope
204
205#ifdef debug
206      G4cout<<"G4QCohChEx::GetMeanFreePath:jI="<<j<<",Zt="<<Z<<",Nt="<<N<<",Mom="<<Momentum
207            <<", XSec="<<CSI/millibarn<<G4endl;
208#endif
209      curIs->second = CSI;
210      susi+=CSI;                            // Make a sum per isotopes
211      SPI->push_back(susi);                 // Remember summed cross-section
212    } // End of temporary initialization of the cross sections in the G4QIsotope singeltone
213    sigma+=Isotopes->GetMeanCrossSection(Z,indEl)*NOfNucPerVolume[i];//SUM(MeanCS*NOfNperV)
214#ifdef debug
215    G4cout<<"G4QCohChEx::GMFP:<S>="<<Isotopes->GetMeanCrossSection(Z,indEl)<<",AddToSigma="
216          <<Isotopes->GetMeanCrossSection(Z,indEl)*NOfNucPerVolume[i]<<G4endl;
217#endif
218    ElProbInMat.push_back(sigma);
219  } // End of LOOP over Elements
220  // Check that cross section is not zero and return the mean free path
221#ifdef debug
222  G4cout<<"G4QCoherentChargeExchange::GetMeanFreePath: MeanFreePath="<<1./sigma<<G4endl;
223#endif
224  if(sigma > 0.) return 1./sigma;                 // Mean path [distance]
225  return DBL_MAX;
226}
227
228G4bool G4QCoherentChargeExchange::IsApplicable(const G4ParticleDefinition& particle) 
229{
230  if      (particle == *(        G4Proton::Proton()        )) return true;
231  else if (particle == *(       G4Neutron::Neutron()       )) return true;
232  //else if (particle == *(     G4MuonMinus::MuonMinus()     )) return true;
233  //else if (particle == *(       G4TauPlus::TauPlus()       )) return true;
234  //else if (particle == *(      G4TauMinus::TauMinus()      )) return true;
235  //else if (particle == *(      G4Electron::Electron()      )) return true;
236  //else if (particle == *(      G4Positron::Positron()      )) return true;
237  //else if (particle == *(         G4Gamma::Gamma()         )) return true;
238  //else if (particle == *(      G4MuonPlus::MuonPlus()      )) return true;
239  //else if (particle == *(G4AntiNeutrinoMu::AntiNeutrinoMu())) return true;
240  //else if (particle == *(   G4NeutrinoMu::NeutrinoMu()   )) return true;
241  //else if (particle == *(     G4PionMinus::PionMinus()     )) return true;
242  //else if (particle == *(      G4PionPlus::PionPlus()      )) return true;
243  //else if (particle == *(      G4KaonPlus::KaonPlus()      )) return true;
244  //else if (particle == *(     G4KaonMinus::KaonMinus()     )) return true;
245  //else if (particle == *(  G4KaonZeroLong::KaonZeroLong()  )) return true;
246  //else if (particle == *( G4KaonZeroShort::KaonZeroShort() )) return true;
247  //else if (particle == *(        G4Lambda::Lambda()        )) return true;
248  //else if (particle == *(     G4SigmaPlus::SigmaPlus()     )) return true;
249  //else if (particle == *(    G4SigmaMinus::SigmaMinus()    )) return true;
250  //else if (particle == *(     G4SigmaZero::SigmaZero()     )) return true;
251  //else if (particle == *(       G4XiMinus::XiMinus()       )) return true;
252  //else if (particle == *(        G4XiZero::XiZero()        )) return true;
253  //else if (particle == *(    G4OmegaMinus::OmegaMinus()    )) return true;
254  //else if (particle == *(   G4AntiNeutron::AntiNeutron()   )) return true;
255  //else if (particle == *(    G4AntiProton::AntiProton()    )) return true;
256#ifdef debug
257  G4cout<<"***>>G4QCoherChargeExch::IsApplicable: PDG="<<particle.GetPDGEncoding()<<G4endl;
258#endif
259  return false;
260}
261
262G4VParticleChange* G4QCoherentChargeExchange::PostStepDoIt(const G4Track& track,
263                                                           const G4Step& step)
264{
265  static const G4double mProt= G4QPDGCode(2212).GetMass();     // CHIPS pMass in MeV
266  static const G4double mNeut= G4QPDGCode(2212).GetMass();     // CHIPS pMass in MeV
267  //
268  //-------------------------------------------------------------------------------------
269  static G4bool CWinit = true;                       // CHIPS Warld needs to be initted
270  if(CWinit)
271  {
272    CWinit=false;
273    G4QCHIPSWorld::Get()->GetParticles(nPartCWorld); // Create CHIPS World (234 part.max)
274  }
275  //-------------------------------------------------------------------------------------
276  const G4DynamicParticle* projHadron = track.GetDynamicParticle();
277  const G4ParticleDefinition* particle=projHadron->GetDefinition();
278#ifdef debug
279  G4cout<<"G4QCohChargeExchange::PostStepDoIt: Before the GetMeanFreePath is called In4M="
280        <<projHadron->Get4Momentum()<<" of PDG="<<particle->GetPDGEncoding()<<", Type="
281        <<particle->GetParticleType()<<", Subtp="<<particle->GetParticleSubType()<<G4endl;
282#endif
283  G4ForceCondition cond=NotForced;
284  GetMeanFreePath(track, -27., &cond);                  // @@ ?? jus to update parameters?
285#ifdef debug
286  G4cout<<"G4QCohChargeExchange::PostStepDoIt: After GetMeanFreePath is called"<<G4endl;
287#endif
288  G4LorentzVector proj4M=(projHadron->Get4Momentum())/MeV; // Convert to MeV!
289  G4double momentum = projHadron->GetTotalMomentum()/MeV; // 3-momentum of the Proj in MeV
290  G4double Momentum = proj4M.rho();                   // @@ Just for the test purposes
291  if(std::fabs(Momentum-momentum)>.000001)
292       G4cerr<<"*Warning*G4QCohChEx::PostStepDoIt:P(IU)="<<Momentum<<"#"<<momentum<<G4endl;
293#ifdef pdebug
294  G4cout<<"G4QCoherentChargeExchange::PostStepDoIt: pP(IU)="<<Momentum<<"="<<momentum
295        <<",pM="<<pM<<",proj4M="<<proj4M<<proj4M.m()<<G4endl;
296#endif
297  if (!IsApplicable(*particle))  // Check applicability
298  {
299    G4cerr<<"G4QCoherentChargeExchange::PostStepDoIt: Only NA is implemented."<<G4endl;
300    return 0;
301  }
302  const G4Material* material = track.GetMaterial();      // Get the current material
303  G4int Z=0;
304  const G4ElementVector* theElementVector = material->GetElementVector();
305  G4int nE=material->GetNumberOfElements();
306#ifdef debug
307  G4cout<<"G4QCohChargeExchange::PostStepDoIt: "<<nE<<" elements in the material."<<G4endl;
308#endif
309  G4int projPDG=0;                           // PDG Code prototype for the captured hadron
310  // Not all these particles are implemented yet (see Is Applicable)
311  if      (particle ==          G4Proton::Proton()         ) projPDG= 2212;
312  else if (particle ==         G4Neutron::Neutron()        ) projPDG= 2112;
313  //else if (particle ==       G4PionMinus::PionMinus()      ) projPDG= -211;
314  //else if (particle ==        G4PionPlus::PionPlus()       ) projPDG=  211;
315  //else if (particle ==        G4KaonPlus::KaonPlus()       ) projPDG= 2112;
316  //else if (particle ==       G4KaonMinus::KaonMinus()      ) projPDG= -321;
317  //else if (particle ==    G4KaonZeroLong::KaonZeroLong()   ) projPDG=  130;
318  //else if (particle ==   G4KaonZeroShort::KaonZeroShort()  ) projPDG=  310;
319  //else if (particle ==        G4MuonPlus::MuonPlus()       ) projPDG=  -13;
320  //else if (particle ==       G4MuonMinus::MuonMinus()      ) projPDG=   13;
321  //else if (particle ==      G4NeutrinoMu::NeutrinoMu()     ) projPDG=   14;
322  //else if (particle ==  G4AntiNeutrinoMu::AntiNeutrinoMu() ) projPDG=  -14;
323  //else if (particle ==        G4Electron::Electron()       ) projPDG=   11;
324  //else if (particle ==        G4Positron::Positron()       ) projPDG=  -11;
325  //else if (particle ==       G4NeutrinoE::NeutrinoE()      ) projPDG=   12;
326  //else if (particle ==   G4AntiNeutrinoE::AntiNeutrinoE()  ) projPDG=  -12;
327  //else if (particle ==           G4Gamma::Gamma()          ) projPDG=   22;
328  //else if (particle ==         G4TauPlus::TauPlus()        ) projPDG=  -15;
329  //else if (particle ==        G4TauMinus::TauMinus()       ) projPDG=   15;
330  //else if (particle ==     G4NeutrinoTau::NeutrinoTau()    ) projPDG=   16;
331  //else if (particle == G4AntiNeutrinoTau::AntiNeutrinoTau()) projPDG=  -16;
332  //else if (particle ==          G4Lambda::Lambda()         ) projPDG= 3122;
333  //else if (particle ==       G4SigmaPlus::SigmaPlus()      ) projPDG= 3222;
334  //else if (particle ==      G4SigmaMinus::SigmaMinus()     ) projPDG= 3112;
335  //else if (particle ==       G4SigmaZero::SigmaZero()      ) projPDG= 3212;
336  //else if (particle ==         G4XiMinus::XiMinus()        ) projPDG= 3312;
337  //else if (particle ==          G4XiZero::XiZero()         ) projPDG= 3322;
338  //else if (particle ==      G4OmegaMinus::OmegaMinus()     ) projPDG= 3334;
339  //else if (particle ==     G4AntiNeutron::AntiNeutron()    ) projPDG=-2112;
340  //else if (particle ==      G4AntiProton::AntiProton()     ) projPDG=-2212;
341#ifdef debug
342  G4int prPDG=particle->GetPDGEncoding();
343  G4cout<<"G4QCohChrgExchange::PostStepDoIt: projPDG="<<projPDG<<", stPDG="<<prPDG<<G4endl;
344#endif
345  if(!projPDG)
346  {
347    G4cerr<<"*Warning*G4QCoherentChargeExchange::PostStepDoIt:UndefinedProjHadron"<<G4endl;
348    return 0;
349  }
350  //G4double pM2=proj4M.m2();        // in MeV^2
351  //G4double pM=std::sqrt(pM2);      // in MeV
352  G4double pM=mNeut;
353  G4int    fPDG=2112;
354  if(projPDG==2112)
355  {
356    pM=mProt;
357    fPDG=2212;
358  }
359  G4double pM2=pM*pM;
360  // Element treatment
361  G4int EPIM=ElProbInMat.size();
362#ifdef debug
363  G4cout<<"G4QCohChEx::PostStDoIt:m="<<EPIM<<",n="<<nE<<",T="<<ElProbInMat[EPIM-1]<<G4endl;
364#endif
365  G4int i=0;
366  if(EPIM>1)
367  {
368    G4double rnd = ElProbInMat[EPIM-1]*G4UniformRand();
369    for(i=0; i<nE; ++i)
370    {
371#ifdef debug
372      G4cout<<"G4QCohChEx::PostStepDoIt:EPM["<<i<<"]="<<ElProbInMat[i]<<",r="<<rnd<<G4endl;
373#endif
374      if (rnd<ElProbInMat[i]) break;
375    }
376    if(i>=nE) i=nE-1;                        // Top limit for the Element
377  }
378  G4Element* pElement=(*theElementVector)[i];
379  Z=static_cast<G4int>(pElement->GetZ());
380#ifdef debug
381    G4cout<<"G4QCoherentChargeExchange::PostStepDoIt: i="<<i<<", Z(element)="<<Z<<G4endl;
382#endif
383  if(Z<=0)
384  {
385    G4cerr<<"-Warning-G4QCoherentChargeExchange::PostStepDoIt: Element with Z="<<Z<<G4endl;
386    if(Z<0) return 0;
387  }
388  std::vector<G4double>* SPI = IsoProbInEl[i];// Vector of summedProbabilities for isotopes
389  std::vector<G4int>* IsN = ElIsoN[i];     // Vector of "#of neutrons" in the isotope El[i]
390  G4int nofIsot=SPI->size();               // #of isotopes in the element i
391#ifdef debug
392  G4cout<<"G4QCohChargeExchange::PosStDoIt:nI="<<nofIsot<<",T="<<(*SPI)[nofIsot-1]<<G4endl;
393#endif
394  G4int j=0;
395  if(nofIsot>1)
396  {
397    G4double rndI=(*SPI)[nofIsot-1]*G4UniformRand(); // Randomize the isotop of the Element
398    for(j=0; j<nofIsot; ++j)
399    {
400#ifdef debug
401      G4cout<<"G4QCohChargEx::PostStepDoIt: SP["<<j<<"]="<<(*SPI)[j]<<", r="<<rndI<<G4endl;
402#endif
403      if(rndI < (*SPI)[j]) break;
404    }
405    if(j>=nofIsot) j=nofIsot-1;            // Top limit for the isotope
406  }
407  G4int N =(*IsN)[j]; ;                    // Randomized number of neutrons
408#ifdef debug
409  G4cout<<"G4QCohChargeEx::PostStepDoIt: j="<<i<<", N(isotope)="<<N<<", MeV="<<MeV<<G4endl;
410#endif
411  if(N<0)
412  {
413    G4cerr<<"*Warning*G4QCohChEx::PostStepDoIt: Isotope with Z="<<Z<<", 0>N="<<N<<G4endl;
414    return 0;
415  }
416  nOfNeutrons=N;                           // Remember it for the energy-momentum check
417#ifdef debug
418  G4cout<<"G4QCohChargeExchange::PostStepDoIt: N="<<N<<" for element with Z="<<Z<<G4endl;
419#endif
420  if(N<0)
421  {
422    G4cerr<<"*Warning*G4QCoherentChargeExchange::PostStepDoIt:Element with N="<<N<< G4endl;
423    return 0;
424  }
425  aParticleChange.Initialize(track);
426#ifdef debug
427  G4cout<<"G4QCoherentChargeExchange::PostStepDoIt: track is initialized"<<G4endl;
428#endif
429  G4double      weight    = track.GetWeight();
430  G4double      localtime = track.GetGlobalTime();
431  G4ThreeVector position  = track.GetPosition();
432#ifdef debug
433  G4cout<<"G4QCoherentChargeExchange::PostStepDoIt: before Touchable extraction"<<G4endl;
434#endif
435  G4TouchableHandle trTouchable = track.GetTouchableHandle();
436#ifdef debug
437  G4cout<<"G4QCoherentChargeExchange::PostStepDoIt: Touchable is extracted"<<G4endl;
438#endif
439  //
440  G4int targPDG=90000000+Z*1000+N;         // CHIPS PDG Code of the target nucleus
441  if(projPDG==2212) targPDG+=999;          // convert to final nucleus code
442  else if(projPDG==2112) targPDG-=999;
443  G4QPDGCode targQPDG(targPDG);            // @@ use G4Ion and get rid of CHIPS World
444  G4double tM=targQPDG.GetMass();          // CHIPS final nucleus mass in MeV
445  G4double kinEnergy= projHadron->GetKineticEnergy()*MeV; // Kin energy in MeV (Is *MeV n?)
446  G4ParticleMomentum dir = projHadron->GetMomentumDirection();// It is a unit three-vector
447  G4LorentzVector tot4M=proj4M+G4LorentzVector(0.,0.,0.,tM); // Total 4-mom of the reaction
448#ifdef debug
449  G4cout<<"G4QCohChrgEx::PostStepDoIt: tM="<<tM<<", p4M="<<proj4M<<", t4M="<<tot4M<<G4endl;
450#endif
451  EnMomConservation=tot4M;                 // Total 4-mom of reaction for E/M conservation
452  // @@ Probably this is not necessary any more
453#ifdef debug
454  G4cout<<"G4QCCX::PSDI:false, P="<<Momentum<<",Z="<<Z<<",N="<<N<<",PDG="<<projPDG<<G4endl;
455#endif
456  G4double xSec=CalculateXSt(false, true, Momentum, Z, N, projPDG); // Recalc. CrossSection
457#ifdef debug
458  G4cout<<"G4QCoChEx::PSDI:PDG="<<projPDG<<",P="<<Momentum<<",CS="<<xSec/millibarn<<G4endl;
459#endif
460#ifdef nandebug
461  if(xSec>0. || xSec<0. || xSec==0);
462  else  G4cout<<"*Warning*G4QCohChargeExchange::PSDI: xSec="<<xSec/millibarn<<G4endl;
463#endif
464  // @@ check a possibility to separate p, n, or alpha (!)
465  if(xSec <= 0.) // The cross-section iz 0 -> Do Nothing
466  {
467#ifdef pdebug
468    G4cerr<<"*Warning*G4QCoherentChargeExchange::PSDoIt:*Zero cross-section* PDG="<<projPDG
469          <<",tPDG="<<targPDG<<",P="<<Momentum<<G4endl;
470#endif
471    //Do Nothing Action insead of the reaction
472    aParticleChange.ProposeEnergy(kinEnergy);
473    aParticleChange.ProposeLocalEnergyDeposit(0.);
474    aParticleChange.ProposeMomentumDirection(dir) ;
475    return G4VDiscreteProcess::PostStepDoIt(track,step);
476  }
477  G4double mint=CalculateXSt(false, false, Momentum, Z, N, projPDG);// randomize t in MeV^2
478#ifdef pdebug
479  G4cout<<"G4QCohChEx::PoStDoIt:pPDG="<<projPDG<<",tPDG="<<targPDG<<",P="<<Momentum<<",CS="
480        <<xSec<<",-t="<<mint<<G4endl;
481#endif
482#ifdef nandebug
483  if(mint>-.0000001);
484  else  G4cout<<"*Warning*G4QCoherentChargeExchange::PostStDoIt:-t="<<mint<<G4endl;
485#endif
486  G4double maxt=CalculateXSt(true, false, Momentum, Z, N, projPDG);
487  if(maxt<=0.) maxt=1.e-27;
488  G4double cost=1.-mint/maxt; // cos(theta) in CMS (@@ we neglect mass diff for ChEx)
489  //
490#ifdef ppdebug
491  G4cout<<"G4QCoherentChargeExchange::PoStDoIt:t="<<mint<<",dpcm2="<<maxt
492        <<",Ek="<<kinEnergy<<",tM="<<tM<<",pM="<<pM<<",cost="<<cost<<G4endl;
493#endif
494  if(cost>1. || cost<-1. || !(cost>-1. || cost<1.))
495  {
496    if(cost>1.000001 || cost<-1.000001 || !(cost>-1. || cost<1.))
497    {
498      G4double tM2=tM*tM;                         // Squared target mass
499      G4double pEn=pM+kinEnergy;                  // tot projectile Energy in MeV
500      G4double sM=(tM+tM)*pEn+tM2+pM2;            // Mondelstam s
501      G4double twop2cm=(tM2+tM2)*(pEn*pEn-pM2)/sM;// Max_t/2 (2*p^2_cm)
502      G4cout<<"*Warning*G4QCohChEx::PostStepDoIt:cos="<<cost<<",t="<<mint<<",T="<<kinEnergy
503            <<",tM="<<tM<<",tmax="<<2*kinEnergy*tM<<",p="<<projPDG<<",t="<<targPDG<<G4endl;
504      G4cout<<"..G4QCohChEx::PoStDoI: dpcm2="<<twop2cm<<"="<<maxt<<G4endl;
505    }
506    if     (cost>1.)  cost=1.;
507    else if(cost<-1.) cost=-1.;
508  }
509  G4LorentzVector reco4M=G4LorentzVector(0.,0.,0.,tM);      // 4mom of the recoil target
510  G4LorentzVector scat4M=G4LorentzVector(0.,0.,0.,pM);      // 4mom of the recoil target
511  G4LorentzVector dir4M=tot4M-G4LorentzVector(0.,0.,0.,(tot4M.e()-tM-pM)*.01);
512  if(!G4QHadron(tot4M).RelDecayIn2(scat4M, reco4M, dir4M, cost, cost))
513  {
514    G4cerr<<"G4QCohChEx::PSDI:t4M="<<tot4M<<",pM="<<pM<<",tM="<<tM<<",cost="<<cost<<G4endl;
515    //G4Exception("G4QCoherentChargeExchange::PostStepDoIt:","009",FatalException,"Decay");
516  }
517#ifdef debug
518  G4cout<<"G4QCohChEx::PoStDoIt:s4M="<<scat4M<<"+r4M="<<reco4M<<"="<<scat4M+reco4M<<G4endl;
519  G4cout<<"G4QCohChEx::PoStDoIt: scatE="<<scat4M.e()-pM<<", recoE="<<reco4M.e()-tM<<",d4M="
520        <<tot4M-scat4M-reco4M<<G4endl;
521#endif
522  // Kill scattered hadron
523  aParticleChange.ProposeTrackStatus(fStopAndKill);
524  // Definition of the scattered nucleon
525  G4DynamicParticle* theSec = new G4DynamicParticle; // A secondary for the recoil hadron
526  G4ParticleDefinition* theDefinition=G4Proton::Proton();
527  if(projPDG==2212) theDefinition=G4Neutron::Neutron();
528  theSec->SetDefinition(theDefinition);
529  EnMomConservation-=scat4M;
530  theSec->Set4Momentum(scat4M);
531  G4Track* aNewTrack = new G4Track(theSec, localtime, position );
532  aNewTrack->SetWeight(weight);                                   //    weighted
533  aNewTrack->SetTouchableHandle(trTouchable);
534  aParticleChange.AddSecondary( aNewTrack );
535  // Filling the recoil nucleus
536  theSec = new G4DynamicParticle; // A secondary for the recoil hadron
537  G4int aA = Z+N;
538#ifdef pdebug
539  G4cout<<"G4QCoherentChargeExchange::PostStepDoIt: Ion Z="<<Z<<", A="<<aA<<G4endl;
540#endif
541  theDefinition=G4ParticleTable::GetParticleTable()->FindIon(Z,aA,0,Z);
542  if(!theDefinition)G4cout<<"*Warning*G4QCohChEx::PostStepDoIt:drop PDG="<<targPDG<<G4endl;
543#ifdef pdebug
544  G4cout<<"G4QCohChEx::PostStepDoIt:RecoilName="<<theDefinition->GetParticleName()<<G4endl;
545#endif
546  theSec->SetDefinition(theDefinition);
547  EnMomConservation-=reco4M;
548#ifdef tdebug
549  G4cout<<"G4QCohChEx::PostSDoIt:"<<targPDG<<reco4M<<reco4M.m()<<EnMomConservation<<G4endl;
550#endif
551  theSec->Set4Momentum(reco4M);
552#ifdef debug
553  G4ThreeVector curD=theSec->GetMomentumDirection();
554  G4double curM=theSec->GetMass();
555  G4double curE=theSec->GetKineticEnergy()+curM;
556  G4cout<<"G4QCohChEx::PostStpDoIt:p="<<curD<<curD.mag()<<",e="<<curE<<",m="<<curM<<G4endl;
557#endif
558  // Make a recoil nucleus
559  aNewTrack = new G4Track(theSec, localtime, position );
560  aNewTrack->SetWeight(weight);                                   //    weighted
561  aNewTrack->SetTouchableHandle(trTouchable);
562  aParticleChange.AddSecondary( aNewTrack );
563#ifdef debug
564    G4cout<<"G4QCoherentChargeExchange::PostStepDoIt:*** PostStepDoIt is done ***"<<G4endl;
565#endif
566  return G4VDiscreteProcess::PostStepDoIt(track, step);
567}
568
569G4double G4QCoherentChargeExchange::CalculateXSt(G4bool oxs, G4bool xst, G4double p,
570                                                 G4int Z, G4int N, G4int pPDG) 
571{
572  static G4bool init=false;
573  static G4bool first=true;
574  static G4VQCrossSection* CSmanager;
575  G4QuasiFreeRatios* qfMan=G4QuasiFreeRatios::GetPointer();
576  if(first)                              // Connection with a singletone
577  {
578    CSmanager=G4QElasticCrossSection::GetPointer();
579    first=false;
580  }
581  G4double res=0.;
582  if(oxs && xst)                         // Only the Cross-Section can be returened
583  {
584    res=CSmanager->GetCrossSection(true, p, Z, N, pPDG); // XS for isotope
585    res*=qfMan->ChExElCoef(p*MeV, Z, N, pPDG);
586  }
587  else if(!oxs && xst)                   // Calculate Cross-Section & prepare differential
588  {
589    res=CSmanager->GetCrossSection(false, p, Z, N, pPDG);// XS for isotope + init t-distr.
590    res*=qfMan->ChExElCoef(p*MeV, Z, N, pPDG);
591    // The XS for the nucleus must be calculated the last
592    init=true;
593  }
594  else if(init)                          // Return t-value for scattering (=G4QElastic)
595  {
596    if(oxs) res=CSmanager->GetHMaxT();   // Calculate the max_t value
597    else res=CSmanager->GetExchangeT(Z, N, pPDG); // fanctionally randomized -t in MeV^2
598  }
599  else G4cout<<"*Warning*G4QCohChrgExchange::CalculateXSt: NotInitiatedScattering"<<G4endl;
600  return res;
601}
Note: See TracBrowser for help on using the repository browser.