source: trunk/source/processes/hadronic/models/chiral_inv_phase_space/interface/src/G4QElastic.cc @ 968

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

update processes

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