Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4QElastic.cc
Go to the documentation of this file.
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$
27 //
28 // ---------------- G4QElastic class -----------------
29 // by Mikhail Kossov, December 2003.
30 // G4QElastic class of the CHIPS Simulation Branch in GEANT4
31 // --------------------------------------------------------------------
32 // Short description: A process for hadron-nucleus elastic scattering.
33 // --------------------------------------------------------------------
34 
35 //#define debug
36 //#define pdebug
37 //#define tdebug
38 //#define nandebug
39 //#define ppdebug
40 
41 #include "G4QElastic.hh"
42 #include "G4SystemOfUnits.hh"
43 #include "G4HadronicDeprecate.hh"
44 
45 
46 // Initialization of static vectors
47 //G4int G4QElastic::nPartCWorld=152; // The#of particles initialized in CHIPS World
48 //G4int G4QElastic::nPartCWorld=122; // The#of particles initialized in CHIPS World
49 G4int G4QElastic::nPartCWorld=85; // The#of particles initialized in CHIPS World Reduced
50 std::vector<G4int> G4QElastic::ElementZ; // Z of the element(i) in theLastCalc
51 std::vector<G4double> G4QElastic::ElProbInMat; // SumProbabilityElements in Material
52 std::vector<std::vector<G4int>*> G4QElastic::ElIsoN; // N of isotope(j) of Element(i)
53 std::vector<std::vector<G4double>*>G4QElastic::IsoProbInEl;//SumProbabIsotopes inElementI
54 
55 // Constructor
56 G4QElastic::G4QElastic(const G4String& processName):
57  G4VDiscreteProcess(processName, fHadronic)
58 {
59  G4HadronicDeprecate("G4QElastic");
60 
61 #ifdef debug
62  G4cout<<"G4QElastic::Constructor is called processName="<<processName<<G4endl;
63 #endif
64  if (verboseLevel>0) G4cout << GetProcessName() << " process is created "<< G4endl;
65  //G4QCHIPSWorld::Get()->GetParticles(nPartCWorld); // Create CHIPS World (234 part. max)
66 }
67 
68 // Destructor
70 
71 
73 
75 
76 // output of the function must be in units of length! L=1/sig_V,sig_V=SUM(n(j,i)*sig(j,i)),
77 // where n(i,j) is a number of nuclei of the isotop j of the element i in V=1(lengtUnit^3)
78 // ********** All CHIPS cross sections are calculated in the surface units ************
80 {
81  *Fc = NotForced;
82  const G4DynamicParticle* incidentParticle = aTrack.GetDynamicParticle();
83  G4ParticleDefinition* incidentParticleDefinition=incidentParticle->GetDefinition();
84  if( !IsApplicable(*incidentParticleDefinition))
85  G4cout<<"*W*G4QElastic::GetMeanFreePath: is called for notImplementedParticle"<<G4endl;
86  // Calculate the mean Cross Section for the set of Elements(*Isotopes) in the Material
87  G4double Momentum = incidentParticle->GetTotalMomentum(); // 3-momentum of the Particle
88 #ifdef debug
89  G4double KinEn = incidentParticle->GetKineticEnergy();
90  G4cout<<"G4QElastic::GetMeanFreePath: kinE="<<KinEn<<",Mom="<<Momentum<<G4endl; // Result
91 #endif
92  const G4Material* material = aTrack.GetMaterial(); // Get the current material
93  const G4double* NOfNucPerVolume = material->GetVecNbOfAtomsPerVolume();
94  const G4ElementVector* theElementVector = material->GetElementVector();
95  G4int nE=material->GetNumberOfElements();
96  G4int pPDG=0;
97  if (incidentParticleDefinition == G4Proton::Proton() ) pPDG=2212;
98  else if(incidentParticleDefinition == G4Neutron::Neutron() ) pPDG=2112;
99  else if(incidentParticleDefinition == G4PionPlus::PionPlus() ) pPDG= 211;
100  else if(incidentParticleDefinition == G4PionMinus::PionMinus() ) pPDG=-211;
101  else if(incidentParticleDefinition == G4KaonPlus::KaonPlus() ) pPDG= 321;
102  else if(incidentParticleDefinition == G4KaonMinus::KaonMinus() ) pPDG=-321;
103  else if(incidentParticleDefinition == G4KaonZeroLong::KaonZeroLong() ) pPDG= 130;
104  else if(incidentParticleDefinition == G4KaonZeroShort::KaonZeroShort() ) pPDG= 310;
105  else if(incidentParticleDefinition == G4Lambda::Lambda() ) pPDG= 3122;
106  else if(incidentParticleDefinition == G4SigmaPlus::SigmaPlus() ) pPDG= 3222;
107  else if(incidentParticleDefinition == G4SigmaMinus::SigmaMinus() ) pPDG= 3112;
108  else if(incidentParticleDefinition == G4SigmaZero::SigmaZero() ) pPDG= 3212;
109  else if(incidentParticleDefinition == G4XiMinus::XiMinus() ) pPDG= 3312;
110  else if(incidentParticleDefinition == G4XiZero::XiZero() ) pPDG= 3322;
111  else if(incidentParticleDefinition == G4OmegaMinus::OmegaMinus() ) pPDG= 3334;
112  else if(incidentParticleDefinition == G4AntiNeutron::AntiNeutron() ) pPDG=-2112;
113  else if(incidentParticleDefinition == G4AntiProton::AntiProton() ) pPDG=-2212;
114  else if(incidentParticleDefinition == G4AntiLambda::AntiLambda() ) pPDG=-3122;
115  else if(incidentParticleDefinition == G4AntiSigmaPlus::AntiSigmaPlus() ) pPDG=-3222;
116  else if(incidentParticleDefinition == G4AntiSigmaMinus::AntiSigmaMinus() ) pPDG=-3112;
117  else if(incidentParticleDefinition == G4AntiSigmaZero::AntiSigmaZero() ) pPDG=-3212;
118  else if(incidentParticleDefinition == G4AntiXiMinus::AntiXiMinus() ) pPDG=-3312;
119  else if(incidentParticleDefinition == G4AntiXiZero::AntiXiZero() ) pPDG=-3322;
120  else if(incidentParticleDefinition == G4AntiOmegaMinus::AntiOmegaMinus() ) pPDG=-3334;
121  else G4cout<<"Warning-G4QElastic::GetMFP:"<<incidentParticleDefinition->GetParticleName()
122  <<" isn't implemented (only SU(3) particles)"<<G4endl;
123 #ifdef pdebug
124  G4cout<<"G4QElastic::GetMeanFreePath:"<<nE<<" Elem's in theMaterial,proj="<<pPDG<<G4endl;
125 #endif
126  G4VQCrossSection* CSmanager=0;
127  G4VQCrossSection* CSmanager2=0;
128  if (pPDG==2212) CSmanager=G4QProtonElasticCrossSection::GetPointer();
129  else if(pPDG==2112) CSmanager=G4QNeutronElasticCrossSection::GetPointer();
130  else if(pPDG== 211) CSmanager=G4QPionPlusElasticCrossSection::GetPointer();
131  else if(pPDG==-211) CSmanager=G4QPionMinusElasticCrossSection::GetPointer();
132  else if(pPDG== 321) CSmanager=G4QKaonPlusElasticCrossSection::GetPointer();
133  else if(pPDG==-321) CSmanager=G4QKaonMinusElasticCrossSection::GetPointer();
134  else if(pPDG== 130 || pPDG== 310)
135  {
138  }
139  else if(pPDG==3222) CSmanager=G4QHyperonPlusElasticCrossSection::GetPointer();
140  else if(pPDG>3110 && pPDG<3335) CSmanager=G4QHyperonElasticCrossSection::GetPointer();
141  else if(pPDG>-3335&&pPDG<-1110) CSmanager=G4QAntiBaryonElasticCrossSection::GetPointer();
142  else G4cout<<"*Warning*G4QElastic::GetMeanFreePath: wrong PDG="<<pPDG<<G4endl;
143  G4QIsotope* Isotopes = G4QIsotope::Get(); // Pointer to the G4QIsotopes singleton
144  G4double sigma=0.; // Sums over elements for the material
145  G4int IPIE=IsoProbInEl.size(); // How many old elements?
146  if(IPIE) for(G4int ip=0; ip<IPIE; ++ip) // Clean up the SumProb's of Isotopes (SPI)
147  {
148  std::vector<G4double>* SPI=IsoProbInEl[ip]; // Pointer to the SPI vector
149  SPI->clear();
150  delete SPI;
151  std::vector<G4int>* IsN=ElIsoN[ip]; // Pointer to the N vector
152  IsN->clear();
153  delete IsN;
154  }
155  ElProbInMat.clear(); // Clean up the SumProb's of Elements (SPE)
156  ElementZ.clear(); // Clear the body vector for Z of Elements
157  IsoProbInEl.clear(); // Clear the body vector for SPI
158  ElIsoN.clear(); // Clear the body vector for N of Isotopes
159  for(G4int i=0; i<nE; ++i)
160  {
161  G4Element* pElement=(*theElementVector)[i]; // Pointer to the current element
162  G4int Z = static_cast<G4int>(pElement->GetZ()); // Z of the Element
163  ElementZ.push_back(Z); // Remember Z of the Element
164  G4int isoSize=0; // The default for the isoVectorLength is 0
165  G4int indEl=0; // Index of non-natural element or 0(default)
166  G4IsotopeVector* isoVector=pElement->GetIsotopeVector(); // Get the predefined IsoVect
167  if(isoVector) isoSize=isoVector->size();// Get size of the existing isotopeVector
168 #ifdef debug
169  G4cout<<"G4QElastic::GetMeanFreePath: isovectorLength="<<isoSize<<G4endl; // Result
170 #endif
171  if(isoSize) // The Element has non-trivial abundance set
172  {
173  indEl=pElement->GetIndex()+1; // Index of the non-trivial element is an order
174 #ifdef debug
175  G4cout<<"G4QEl::GetMFP: iE="<<indEl<<", def="<<Isotopes->IsDefined(Z,indEl)<<G4endl;
176 #endif
177  if(!Isotopes->IsDefined(Z,indEl)) // This index is not defined for this Z: define
178  {
179  std::vector<std::pair<G4int,G4double>*>* newAbund =
180  new std::vector<std::pair<G4int,G4double>*>;
181  G4double* abuVector=pElement->GetRelativeAbundanceVector();
182  for(G4int j=0; j<isoSize; j++) // Calculation of abundance vector for isotopes
183  {
184  G4int N=pElement->GetIsotope(j)->GetN()-Z; // N means A=N+Z !
185  if(pElement->GetIsotope(j)->GetZ()!=Z)G4cerr<<"G4QElastic::GetMeanFreePath: Z="
186  <<pElement->GetIsotope(j)->GetZ()<<"#"<<Z<<G4endl;
187  G4double abund=abuVector[j];
188  std::pair<G4int,G4double>* pr= new std::pair<G4int,G4double>(N,abund);
189 #ifdef debug
190  G4cout<<"G4QElastic::GetMeanFreePath:pair#="<<j<<",N="<<N<<",ab="<<abund<<G4endl;
191 #endif
192  newAbund->push_back(pr);
193  }
194 #ifdef debug
195  G4cout<<"G4QElastic::GetMeanFreePath: pairVectorLength="<<newAbund->size()<<G4endl;
196 #endif
197  indEl=G4QIsotope::Get()->InitElement(Z,indEl,newAbund); // definition of the newInd
198  for(G4int k=0; k<isoSize; k++) delete (*newAbund)[k]; // Cleaning temporary
199  delete newAbund; // Was "new" in the beginning of the name space
200  }
201  }
202  std::vector<std::pair<G4int,G4double>*>* cs= Isotopes->GetCSVector(Z,indEl);//CSPointer
203  std::vector<G4double>* SPI = new std::vector<G4double>; // Pointer to the SPI vector
204  IsoProbInEl.push_back(SPI);
205  std::vector<G4int>* IsN = new std::vector<G4int>; // Pointer to the N vector
206  ElIsoN.push_back(IsN);
207  G4int nIs=cs->size(); // A#Of Isotopes in the Element
208 #ifdef debug
209  G4cout<<"G4QEl::GMFP:=***=>,#isot="<<nIs<<", Z="<<Z<<", indEl="<<indEl<<G4endl;
210 #endif
211  G4double susi=0.; // sum of CS over isotopes
212  if(nIs) for(G4int j=0; j<nIs; j++) // Calculate CS for eachIsotope of El
213  {
214  std::pair<G4int,G4double>* curIs=(*cs)[j]; // A pointer, which is used twice
215  G4int N=curIs->first; // #of Neuterons in the isotope j of El i
216  IsN->push_back(N); // Remember Min N for the Element
217 #ifdef debug
218  G4cout<<"G4QEl::GMFP:*true*,P="<<Momentum<<",Z="<<Z<<",N="<<N<<",PDG="<<pPDG<<G4endl;
219 #endif
220  G4bool ccsf=true;
221  if(Q==-27.) ccsf=false;
222 #ifdef debug
223  G4cout<<"G4QEl::GMFP: GetCS #1 j="<<j<<G4endl;
224 #endif
225  G4double CSI=0.; // Prototype of CS(j,i) for the isotope
226  if(!CSmanager2) CSI=CSmanager->GetCrossSection(ccsf,Momentum,Z,N,pPDG); // CS(j,i)
227  else CSI=(CSmanager->GetCrossSection(ccsf,Momentum,Z,N,-321)+
228  CSmanager2->GetCrossSection(ccsf,Momentum,Z,N,321) )/2.; // K0
229 #ifdef debug
230  G4cout<<"G4QEl::GMFP: jI="<<j<<", Zt="<<Z<<", Nt="<<N<<", Mom="<<Momentum<<", XSec="
231  <<CSI/millibarn<<G4endl;
232 #endif
233  curIs->second = CSI;
234  susi+=CSI; // Make a sum per isotopes
235  SPI->push_back(susi); // Remember summed cross-section
236  } // End of temporary initialization of the cross sections in the G4QIsotope singeltone
237  sigma+=Isotopes->GetMeanCrossSection(Z,indEl)*NOfNucPerVolume[i];//SUM(MeanCS*NOfNperV)
238 #ifdef debug
239  G4cout<<"G4QEl::GMFP: <S>="<<Isotopes->GetMeanCrossSection(Z,indEl)<<", AddToSigma="
240  <<Isotopes->GetMeanCrossSection(Z,indEl)*NOfNucPerVolume[i]<<G4endl;
241 #endif
242  ElProbInMat.push_back(sigma);
243  } // End of LOOP over Elements
244  // Check that cross section is not zero and return the mean free path
245 #ifdef debug
246  G4cout<<"G4QElastic::GetMeanFreePath: MeanFreePath="<<1./sigma<<G4endl;
247 #endif
248  if(sigma > 0.) return 1./sigma; // Mean path [distance]
249  return DBL_MAX;
250 }
251 
252 
254 {
255  if (particle == *( G4Proton::Proton() )) return true;
256  else if (particle == *( G4Neutron::Neutron() )) return true;
257  else if (particle == *( G4MuonMinus::MuonMinus() )) return true;
258  else if (particle == *( G4TauPlus::TauPlus() )) return true;
259  else if (particle == *( G4TauMinus::TauMinus() )) return true;
260  else if (particle == *( G4Electron::Electron() )) return true;
261  else if (particle == *( G4Positron::Positron() )) return true;
262  else if (particle == *( G4Gamma::Gamma() )) return true;
263  else if (particle == *( G4MuonPlus::MuonPlus() )) return true;
264  else if (particle == *(G4AntiNeutrinoMu::AntiNeutrinoMu())) return true;
265  else if (particle == *( G4NeutrinoMu::NeutrinoMu() )) return true;
266  else if (particle == *( G4PionMinus::PionMinus() )) return true;
267  else if (particle == *( G4PionPlus::PionPlus() )) return true;
268  else if (particle == *( G4KaonPlus::KaonPlus() )) return true;
269  else if (particle == *( G4KaonMinus::KaonMinus() )) return true;
270  else if (particle == *( G4KaonZeroLong::KaonZeroLong() )) return true;
271  else if (particle == *( G4KaonZeroShort::KaonZeroShort() )) return true;
272  else if (particle == *( G4Lambda::Lambda() )) return true;
273  else if (particle == *( G4SigmaPlus::SigmaPlus() )) return true;
274  else if (particle == *( G4SigmaMinus::SigmaMinus() )) return true;
275  else if (particle == *( G4SigmaZero::SigmaZero() )) return true;
276  else if (particle == *( G4XiMinus::XiMinus() )) return true;
277  else if (particle == *( G4XiZero::XiZero() )) return true;
278  else if (particle == *( G4OmegaMinus::OmegaMinus() )) return true;
279  else if (particle == *( G4AntiNeutron::AntiNeutron() )) return true;
280  else if (particle == *( G4AntiProton::AntiProton() )) return true;
281 #ifdef debug
282  G4cout<<"***>>G4QElastic::IsApplicable: PDG="<<particle.GetPDGEncoding()<<G4endl;
283 #endif
284  return false;
285 }
286 
288 {
289  //static const G4double mProt=G4Proton::Proton()->GetPDGMass()*GeV;// proton mass in GeV
290  //static const G4double mProt= G4QPDGCode(2212).GetMass()*.001; // CHIPS m_p in GeV
291  //static const G4double mP2=mProt*mProt; // squared proton mass
292  //
293  //-------------------------------------------------------------------------------------
294  static G4bool CWinit = true; // CHIPS Warld needs to be initted
295  if(CWinit)
296  {
297  CWinit=false;
298  G4QCHIPSWorld::Get()->GetParticles(nPartCWorld); // Create CHIPS World (234 part.max)
299  }
300  //-------------------------------------------------------------------------------------
301  const G4DynamicParticle* projHadron = track.GetDynamicParticle();
302  const G4ParticleDefinition* particle=projHadron->GetDefinition();
303 #ifdef debug
304  G4cout<<"G4QElastic::PostStepDoIt: Before the GetMeanFreePath is called In4M="
305  <<projHadron->Get4Momentum()<<" of PDG="<<particle->GetPDGEncoding()<<", Type="
306  <<particle->GetParticleType()<<", Subtp="<<particle->GetParticleSubType()<<G4endl;
307 #endif
309  GetMeanFreePath(track, -27., &cond); // @@ ?? jus to update parameters?
310 #ifdef debug
311  G4cout<<"G4QElastic::PostStepDoIt: After the GetMeanFreePath is called"<<G4endl;
312 #endif
313  G4LorentzVector proj4M=(projHadron->Get4Momentum())/MeV; // Convert to MeV!
314  G4LorentzVector scat4M=proj4M; // @@ Must be filled (?)
315  G4double momentum = projHadron->GetTotalMomentum()/MeV; // 3-momentum of the Proj in MeV
316  G4double Momentum = proj4M.rho(); // @@ Just for the test purposes
317  if(std::fabs(Momentum-momentum)>.000001)
318  G4cerr<<"*War*G4QElastic::PostStepDoIt:P(IU)="<<Momentum<<"#"<<momentum<<G4endl;
319  G4double pM2=proj4M.m2(); // in MeV^2
320  G4double pM=std::sqrt(pM2); // in MeV
321 #ifdef debug
322  G4cout<<"G4QElastic::PostStepDoIt: pP(IU)="<<Momentum<<"="<<momentum<<",pM="<<pM
323  <<",scat4M="<<scat4M<<scat4M.m()<<G4endl;
324 #endif
325  if (!IsApplicable(*particle)) // Check applicability
326  {
327  G4cerr<<"G4QElastic::PostStepDoIt: Only NA elastic is implemented."<<G4endl;
328  return 0;
329  }
330  const G4Material* material = track.GetMaterial(); // Get the current material
331  const G4ElementVector* theElementVector = material->GetElementVector();
332  G4int nE=material->GetNumberOfElements();
333 #ifdef debug
334  G4cout<<"G4QElastic::PostStepDoIt: "<<nE<<" elements in the material."<<G4endl;
335 #endif
336  G4int projPDG=0; // PDG Code prototype for the captured hadron
337  // Not all these particles are implemented yet (see Is Applicable)
338  if (particle == G4Proton::Proton() ) projPDG= 2212;
339  else if (particle == G4Neutron::Neutron() ) projPDG= 2112;
340  else if (particle == G4PionMinus::PionMinus() ) projPDG= -211;
341  else if (particle == G4PionPlus::PionPlus() ) projPDG= 211;
342  else if (particle == G4KaonPlus::KaonPlus() ) projPDG= 2112;
343  else if (particle == G4KaonMinus::KaonMinus() ) projPDG= -321;
344  else if (particle == G4KaonZeroLong::KaonZeroLong() ) projPDG= 130;
345  else if (particle == G4KaonZeroShort::KaonZeroShort() ) projPDG= 310;
346  //else if (particle == G4MuonPlus::MuonPlus() ) projPDG= -13;
347  //else if (particle == G4MuonMinus::MuonMinus() ) projPDG= 13;
348  //else if (particle == G4NeutrinoMu::NeutrinoMu() ) projPDG= 14;
349  //else if (particle == G4AntiNeutrinoMu::AntiNeutrinoMu() ) projPDG= -14;
350  //else if (particle == G4Electron::Electron() ) projPDG= 11;
351  //else if (particle == G4Positron::Positron() ) projPDG= -11;
352  //else if (particle == G4NeutrinoE::NeutrinoE() ) projPDG= 12;
353  //else if (particle == G4AntiNeutrinoE::AntiNeutrinoE() ) projPDG= -12;
354  //else if (particle == G4Gamma::Gamma() ) projPDG= 22;
355  //else if (particle == G4TauPlus::TauPlus() ) projPDG= -15;
356  //else if (particle == G4TauMinus::TauMinus() ) projPDG= 15;
357  //else if (particle == G4NeutrinoTau::NeutrinoTau() ) projPDG= 16;
358  //else if (particle == G4AntiNeutrinoTau::AntiNeutrinoTau()) projPDG= -16;
359  else if (particle == G4Lambda::Lambda() ) projPDG= 3122;
360  else if (particle == G4SigmaPlus::SigmaPlus() ) projPDG= 3222;
361  else if (particle == G4SigmaMinus::SigmaMinus() ) projPDG= 3112;
362  else if (particle == G4SigmaZero::SigmaZero() ) projPDG= 3212;
363  else if (particle == G4XiMinus::XiMinus() ) projPDG= 3312;
364  else if (particle == G4XiZero::XiZero() ) projPDG= 3322;
365  else if (particle == G4OmegaMinus::OmegaMinus() ) projPDG= 3334;
366  else if (particle == G4AntiNeutron::AntiNeutron() ) projPDG=-2112;
367  else if (particle == G4AntiProton::AntiProton() ) projPDG=-2212;
368  else if (particle == G4AntiLambda::AntiLambda() ) projPDG=-3122;
369  else if (particle == G4AntiSigmaPlus::AntiSigmaPlus() ) projPDG=-3222;
370  else if (particle == G4AntiSigmaMinus::AntiSigmaMinus() ) projPDG=-3112;
371  else if (particle == G4AntiSigmaZero::AntiSigmaZero() ) projPDG=-3212;
372  else if (particle == G4AntiXiMinus::AntiXiMinus() ) projPDG=-3312;
373  else if (particle == G4AntiXiZero::AntiXiZero() ) projPDG=-3322;
374  else if (particle == G4AntiOmegaMinus::AntiOmegaMinus() ) projPDG=-3334;
375 #ifdef pdebug
376  G4int prPDG=particle->GetPDGEncoding();
377  G4cout<<"G4QElastic::PostStepDoIt: projPDG="<<projPDG<<", stPDG="<<prPDG<<G4endl;
378 #endif
379  if(!projPDG)
380  {
381  G4cerr<<"*Warning*G4QElastic::PostStepDoIt: Undefined interacting hadron"<<G4endl;
382  return 0;
383  }
384  G4int EPIM=ElProbInMat.size();
385 #ifdef debug
386  G4cout<<"G4QElastic::PostStDoIt:m="<<EPIM<<",n="<<nE<<",T="<<ElProbInMat[EPIM-1]<<G4endl;
387 #endif
388  G4int i=0;
389  if(EPIM>1)
390  {
391  G4double rnd = ElProbInMat[EPIM-1]*G4UniformRand();
392  for(i=0; i<nE; ++i)
393  {
394 #ifdef debug
395  G4cout<<"G4QElastic::PostStepDoIt:EPM["<<i<<"]="<<ElProbInMat[i]<<",r="<<rnd<<G4endl;
396 #endif
397  if (rnd<ElProbInMat[i]) break;
398  }
399  if(i>=nE) i=nE-1; // Top limit for the Element
400  }
401  G4Element* pElement=(*theElementVector)[i];
402  G4int Z=static_cast<G4int>(pElement->GetZ());
403 #ifdef debug
404  G4cout<<"G4QElastic::PostStepDoIt: i="<<i<<", Z(element)="<<Z<<G4endl;
405 #endif
406  if(Z<=0)
407  {
408  G4cerr<<"---Warning---G4QElastic::PostStepDoIt: Element with Z="<<Z<<G4endl;
409  if(Z<0) return 0;
410  }
411  std::vector<G4double>* SPI = IsoProbInEl[i];// Vector of summedProbabilities for isotopes
412  std::vector<G4int>* IsN = ElIsoN[i]; // Vector of "#of neutrons" in the isotope El[i]
413  G4int nofIsot=SPI->size(); // #of isotopes in the element i
414 #ifdef debug
415  G4cout<<"G4QElastic::PosStDoIt: nI="<<nofIsot<<",T="<<(*SPI)[nofIsot-1]<<G4endl;
416 #endif
417  G4int j=0;
418  if(nofIsot>1)
419  {
420  G4double rndI=(*SPI)[nofIsot-1]*G4UniformRand(); // Randomize the isotop of the Element
421  for(j=0; j<nofIsot; ++j)
422  {
423 #ifdef debug
424  G4cout<<"G4QElastic::PostStepDoIt: SP["<<j<<"]="<<(*SPI)[j]<<", r="<<rndI<<G4endl;
425 #endif
426  if(rndI < (*SPI)[j]) break;
427  }
428  if(j>=nofIsot) j=nofIsot-1; // Top limit for the isotope
429  }
430  G4int N =(*IsN)[j]; ; // Randomized number of neutrons
431 #ifdef debug
432  G4cout<<"G4QElastic::PostStepDoIt: j="<<i<<", N(isotope)="<<N<<", MeV="<<MeV<<G4endl;
433 #endif
434  if(N<0)
435  {
436  G4cerr<<"*Warning*G4QElastic::PostStepDoIt: Isotope with Z="<<Z<<", 0>N="<<N<<G4endl;
437  return 0;
438  }
439  nOfNeutrons=N; // Remember it for the energy-momentum check
440 #ifdef debug
441  G4cout<<"G4QElastic::PostStepDoIt: N="<<N<<" for element with Z="<<Z<<G4endl;
442 #endif
444 #ifdef debug
445  G4cout<<"G4QElastic::PostStepDoIt: track is initialized"<<G4endl;
446 #endif
447  G4double weight = track.GetWeight();
448  G4double localtime = track.GetGlobalTime();
450 #ifdef debug
451  G4cout<<"G4QElastic::PostStepDoIt: before Touchable extraction"<<G4endl;
452 #endif
453  G4TouchableHandle trTouchable = track.GetTouchableHandle();
454 #ifdef debug
455  G4cout<<"G4QElastic::PostStepDoIt: Touchable is extracted"<<G4endl;
456 #endif
457  //
458  G4int targPDG=90000000+Z*1000+N; // CHIPS PDG Code of the target nucleus
459  G4QPDGCode targQPDG(targPDG); // @@ one can use G4Ion & get rid of CHIPSWorld
460  G4double tM=targQPDG.GetMass(); // CHIPS target mass in MeV
461  G4double kinEnergy= projHadron->GetKineticEnergy()*MeV; // Kin energy in MeV (Is *MeV n?)
462  G4ParticleMomentum dir = projHadron->GetMomentumDirection();// It is a unit three-vector
463  G4LorentzVector tot4M=proj4M+G4LorentzVector(0.,0.,0.,tM); // Total 4-mom of the reaction
464 #ifdef debug
465  G4cout<<"G4QElastic::PostStepDoIt: tM="<<tM<<", p4M="<<proj4M<<", t4M="<<tot4M<<G4endl;
466 #endif
467  EnMomConservation=tot4M; // Total 4-mom of reaction for E/M conservation
468  G4VQCrossSection* CSmanager=0;
469  G4VQCrossSection* CSmanager2=0;
470  if (projPDG==2212) CSmanager=G4QProtonElasticCrossSection::GetPointer();
471  else if(projPDG==2112) CSmanager=G4QNeutronElasticCrossSection::GetPointer();
472  else if(projPDG== 211) CSmanager=G4QPionPlusElasticCrossSection::GetPointer();
473  else if(projPDG==-211) CSmanager=G4QPionMinusElasticCrossSection::GetPointer();
474  else if(projPDG== 321) CSmanager=G4QKaonPlusElasticCrossSection::GetPointer();
475  else if(projPDG==-321) CSmanager=G4QKaonMinusElasticCrossSection::GetPointer();
476  else if(projPDG== 130 || projPDG== 310)
477  {
480  }
481  else if(projPDG==3222) CSmanager=G4QHyperonPlusElasticCrossSection::GetPointer();
482  else if(projPDG>3110&&projPDG<3335)CSmanager=G4QHyperonElasticCrossSection::GetPointer();
483  else if(projPDG>-3335 && projPDG<-1110)
485  else G4cout<<"*Warning*G4QElastic::PostStepDoIt: wrong PDG="<<projPDG<<G4endl;
486 #ifdef debug
487  G4cout<<"G4QElas::PSDI:false,P="<<Momentum<<",Z="<<Z<<",N="<<N<<",PDG="<<projPDG<<G4endl;
488 #endif
489  G4double xSec=0.;
490  if(!CSmanager2) xSec=CSmanager->GetCrossSection(false,Momentum,Z,N,projPDG); //ReactionCS
491  else xSec=(CSmanager->GetCrossSection(false,Momentum,Z,N,-321)+
492  CSmanager2->GetCrossSection(false,Momentum,Z,N,321) )/2.; // K0
493 #ifdef debug
494  G4cout<<"G4QElast::PSDI:pPDG="<<projPDG<<",P="<<Momentum<<",CS="<<xSec/millibarn<<G4endl;
495 #endif
496 #ifdef nandebug
497  if(xSec>0. || xSec<0. || xSec==0);
498  else G4cout<<"*Warning*G4QElastic::PostStDoIt: xSec="<<xSec/millibarn<<G4endl;
499 #endif
500  // @@ check a possibility to separate p, n, or alpha (!)
501  if(xSec <= 0.) // The cross-section iz 0 -> Do Nothing
502  {
503 #ifdef debug
504  G4cerr<<"*Warning*G4QElastic::PostStDoIt: *Zero cross-section* PDG="<<projPDG<<",tPDG="
505  <<targPDG<<",P="<<Momentum<<G4endl;
506 #endif
507  //Do Nothing Action insead of the reaction
508  aParticleChange.ProposeEnergy(kinEnergy);
511  return G4VDiscreteProcess::PostStepDoIt(track,step);
512  }
513  G4double mint=CSmanager->GetExchangeT(Z, N, projPDG); // functional rand -t(MeV^2)
514 #ifdef debug
515  G4cout<<"G4QElast::PSDI:pPDG="<<projPDG<<",tPDG="<<targPDG<<",P="<<Momentum<<",CS="
516  <<xSec<<",-t="<<mint<<G4endl;
517 #endif
518 #ifdef nandebug
519  if(mint>-.0000001);
520  else G4cout<<"*Warning*G4QElastic::PostStDoIt:-t="<<mint<<G4endl;
521 #endif
522  G4double maxt=CSmanager->GetHMaxT(); // max -t (MeV^2)
523  G4double cost=1.-mint/maxt; // cos(theta) in CMS
524  //
525 #ifdef ppdebug
526  G4cout<<"G4QElastic::PoStDoI:t="<<mint<<", maxt="<<maxt<<",Ek="<<kinEnergy<<",tM="<<tM
527  <<",pM="<<pM<<",cost="<<cost<<G4endl;
528 #endif
529  if(cost>1. || cost<-1. || !(cost>-1. || cost<1.))
530  {
531  if(cost>1.000001 || cost<-1.000001 || !(cost>-1. || cost<1.))
532  {
533  G4double tM2=tM*tM; // Squared target mass
534  G4double pEn=pM+kinEnergy; // tot projectile Energy in MeV
535  G4double sM=(tM+tM)*pEn+tM2+pM2; // Mondelstam s
536  G4double twop2cm=(tM2+tM2)*(pEn*pEn-pM2)/sM;// Max_t/2 (2*p^2_cm)
537  G4cout<<"*Warning*G4QElastic::PostStepDoIt:cos="<<cost<<",t="<<mint<<",T="<<kinEnergy
538  <<",tM="<<tM<<",tmax="<<2*kinEnergy*tM<<",p="<<projPDG<<",t="<<targPDG<<G4endl;
539  G4cout<<"..G4QElastic::PoStDoI: dpcm2="<<twop2cm<<"="<<maxt<<G4endl;
540  }
541  if (cost>1.) cost=1.;
542  else if(cost<-1.) cost=-1.;
543  }
544  G4LorentzVector reco4M=G4LorentzVector(0.,0.,0.,tM); // 4mom of the recoil target
545  G4LorentzVector dir4M=tot4M-G4LorentzVector(0.,0.,0.,(tot4M.e()-tM-pM)*.01);
546  if(!G4QHadron(tot4M).RelDecayIn2(scat4M, reco4M, dir4M, cost, cost))
547  {
548  G4cerr<<"G4QElastic::PSD:t4M="<<tot4M<<",pM="<<pM<<",tM="<<tM<<",cost="<<cost<<G4endl;
549  //G4Exception("G4QElastic::PostStepDoIt:","009",FatalException,"Decay of ElasticComp");
550  }
551 #ifdef debug
552  G4cout<<"G4QElastic::PoStDoIt:s4M="<<scat4M<<"+r4M="<<reco4M<<"="<<scat4M+reco4M<<G4endl;
553  G4cout<<"G4QElastic::PoStDoIt: scatE="<<scat4M.e()-pM<<", recoE="<<reco4M.e()-tM<<",d4M="
554  <<tot4M-scat4M-reco4M<<G4endl;
555 #endif
556  // Update G4VParticleChange for the scattered projectile
557  G4double finE=scat4M.e()-pM; // Final kinetic energy of the scattered proton
558  if(finE>0.0) aParticleChange.ProposeEnergy(finE);
559  else
560  {
561  if(finE<-1.e-8 || !(finE>-1.||finE<1.)) // NAN or negative
562  G4cerr<<"*Warning*G4QElastic::PostStDoIt: Zero or negative scattered E="<<finE
563  <<", s4M="<<scat4M<<", r4M="<<reco4M<<", d4M="<<tot4M-scat4M-reco4M<<G4endl;
564  //G4Exception("G4QElastic::PostStDoIt()","009", FatalException," <0 or nan energy");
567  }
568  G4ThreeVector findir=scat4M.vect()/scat4M.rho(); // Unit vector in new direction
569  aParticleChange.ProposeMomentumDirection(findir); // new direction for the scattered part
570  EnMomConservation-=scat4M; // It must be initialized by (pE+tM,pP)
571  // This is how in general the secondary should be identified
572  G4DynamicParticle* theSec = new G4DynamicParticle; // A secondary for the recoil hadron
573  G4int aA = Z+N;
574 #ifdef debug
575  G4cout<<"G4QElastic::PostStepDoIt: Ion Z="<<Z<<", A="<<aA<<G4endl;
576 #endif
578  ->FindIon(Z,aA,0,Z);
579  if(!theDefinition)G4cout<<"*Warning*G4QElastic::PostStepDoIt:drop PDG="<<targPDG<<G4endl;
580 #ifdef debug
581  G4cout<<"G4QElastic::PostStepDoIt:RecoilName="<<theDefinition->GetParticleName()<<G4endl;
582 #endif
583  theSec->SetDefinition(theDefinition);
584 
585  EnMomConservation-=reco4M;
586 #ifdef tdebug
587  G4cout<<"G4QElastic::PostSDoIt:"<<targPDG<<reco4M<<reco4M.m()<<EnMomConservation<<G4endl;
588 #endif
589  theSec->Set4Momentum(reco4M);
590 #ifdef debug
591  G4ThreeVector curD=theSec->GetMomentumDirection();
592  G4double curM=theSec->GetMass();
593  G4double curE=theSec->GetKineticEnergy()+curM;
594  G4cout<<"G4QElastic::PSDoIt:p="<<curD<<curD.mag()<<",e="<<curE<<",m="<<curM<<G4endl;
595 #endif
596  // Make a recoil nucleus
597  G4Track* aNewTrack = new G4Track(theSec, localtime, position );
598  aNewTrack->SetWeight(weight); // weighted
599  aNewTrack->SetTouchableHandle(trTouchable);
600  aParticleChange.AddSecondary( aNewTrack );
601 #ifdef debug
602  G4cout<<"G4QElastic::PostStepDoIt: **** PostStepDoIt is done ****"<<G4endl;
603 #endif
604  return G4VDiscreteProcess::PostStepDoIt(track, step);
605 }
606 
607 // @@ For future improvement @@
608 G4double G4QElastic::CalculateXSt(G4bool oxs, G4bool xst, G4double p, G4int Z, G4int N,
609  G4int pPDG)
610 {
611  static G4bool init=false;
612  static G4bool first=true;
613  static G4VQCrossSection* PCSmanager;
614  static G4VQCrossSection* NCSmanager;
615  static G4VQCrossSection* PiPCSmanager;
616  static G4VQCrossSection* PiMCSmanager;
617  static G4VQCrossSection* KaPCSmanager;
618  static G4VQCrossSection* KaMCSmanager;
619  static G4VQCrossSection* HypCSmanager;
620  static G4VQCrossSection* HyPCSmanager;
621  static G4VQCrossSection* aBaCSmanager;
622  if(first) // Connection with a singletone
623  {
633  first=false;
634  }
635  G4double res=0.;
636  if(oxs && xst) // Only the Cross-Section can be returened
637  {
638  if (pPDG==2212) res=PCSmanager->GetCrossSection(true, p, Z, N, pPDG);
639  else if(pPDG==2112) res=NCSmanager->GetCrossSection(true, p, Z, N, pPDG);
640  else if(pPDG== 211) res=PiPCSmanager->GetCrossSection(true, p, Z, N, pPDG);
641  else if(pPDG==-211) res=PiMCSmanager->GetCrossSection(true, p, Z, N, pPDG);
642  else if(pPDG== 321) res=KaPCSmanager->GetCrossSection(true, p, Z, N, pPDG);
643  else if(pPDG==-321) res=KaMCSmanager->GetCrossSection(true, p, Z, N, pPDG);
644  else if(pPDG==310||pPDG==130) res=(KaMCSmanager->GetCrossSection(true, p, Z, N, pPDG)+
645  KaPCSmanager->GetCrossSection(true, p, Z, N, pPDG))/2.;
646  else if(pPDG==3222) res=HyPCSmanager->GetCrossSection(true, p, Z, N, pPDG);
647  else if(pPDG>3110 && pPDG<3335) res=HypCSmanager->GetCrossSection(true, p, Z, N, pPDG);
648  else if(pPDG>-3335 && pPDG<-1110) res=aBaCSmanager->GetCrossSection(true,p,Z,N,pPDG);
649  else G4cout<<"*Warning*G4QElastic::CalculateXSt: (o) wrong PDG="<<pPDG<<G4endl;
650  }
651  else if(!oxs && xst) // Calculate CrossSection & prepare differentialCS
652  {
653  if (pPDG==2212) res=PCSmanager->GetCrossSection(false, p, Z, N, pPDG);
654  else if(pPDG==2112) res=NCSmanager->GetCrossSection(false, p, Z, N, pPDG);
655  else if(pPDG== 211) res=PiPCSmanager->GetCrossSection(false, p, Z, N, pPDG);
656  else if(pPDG==-211) res=PiMCSmanager->GetCrossSection(false, p, Z, N, pPDG);
657  else if(pPDG== 321) res=KaPCSmanager->GetCrossSection(false, p, Z, N, pPDG);
658  else if(pPDG==-321) res=KaMCSmanager->GetCrossSection(false, p, Z, N, pPDG);
659  else if(pPDG==310||pPDG==130) res=(KaMCSmanager->GetCrossSection(false, p, Z, N, pPDG)+
660  KaPCSmanager->GetCrossSection(false, p, Z, N, pPDG))/2.;
661  else if(pPDG==3222) res=HyPCSmanager->GetCrossSection(false, p, Z, N, pPDG);
662  else if(pPDG>3110 && pPDG<3335) res=HypCSmanager->GetCrossSection(false, p, Z,N, pPDG);
663  else if(pPDG>-3335 && pPDG<-1110) res=aBaCSmanager->GetCrossSection(false,p,Z,N,pPDG);
664  else G4cout<<"*Warning*G4QElastic::CalculateXSt: (x) wrong PDG="<<pPDG<<G4endl;
665  // The XS for the nucleus must be calculated the last
666  init=true;
667  }
668  else if(init) // Return t-value for scattering (=G4QElastic)
669  {
670  if (pPDG==2212) // ===> Protons
671  {
672  if(oxs) res=PCSmanager->GetHMaxT(); // Calculate the max_t value
673  else res=PCSmanager->GetExchangeT(Z, N, pPDG); // functionally randomized -t in MeV^2
674  }
675  else if(pPDG==2112) // ==> Neutrons
676  {
677  if(oxs) res=NCSmanager->GetHMaxT(); // Calculate the max_t value
678  else res=NCSmanager->GetExchangeT(Z, N, pPDG); // functionally randomized -t in MeV^2
679  }
680  else if(pPDG== 211) // ==> PionPlus
681  {
682  if(oxs) res=PiPCSmanager->GetHMaxT(); // Calculate the max_t value
683  else res=PiPCSmanager->GetExchangeT(Z, N, pPDG);// functionallyRandomized -t in MeV^2
684  }
685  else if(pPDG==-211) // ==> PionMinus
686  {
687  if(oxs) res=PiMCSmanager->GetHMaxT(); // Calculate the max_t value
688  else res=PiMCSmanager->GetExchangeT(Z, N, pPDG);// functionallyRandomized -t in MeV^2
689  }
690  else if(pPDG==321 || pPDG==310 || pPDG==130) // ==> KaonPlus or KaonZero
691  {
692  if(oxs) res=KaPCSmanager->GetHMaxT(); // Calculate the max_t value
693  else res=KaPCSmanager->GetExchangeT(Z, N, pPDG);// functionallyRandomized -t in MeV^2
694  }
695  else if(pPDG==-321) // ==> KaonMinus
696  {
697  if(oxs) res=KaMCSmanager->GetHMaxT(); // Calculate the max_t value
698  else res=KaMCSmanager->GetExchangeT(Z, N, pPDG);// functionallyRandomized -t in MeV^2
699  }
700  else if(pPDG==3222) // ==> SigmaPlus
701  {
702  if(oxs) res=HyPCSmanager->GetHMaxT(); // Calculate the max_t value
703  else res=HyPCSmanager->GetExchangeT(Z, N, pPDG);// functionallyRandomized -t in MeV^2
704  }
705  else if(pPDG<3335 && pPDG>1110) // ==> Rest of Hyperons
706  {
707  if(oxs) res=HypCSmanager->GetHMaxT(); // Calculate the max_t value
708  else res=HypCSmanager->GetExchangeT(Z, N, pPDG);// functionallyRandomized -t in MeV^2
709  }
710  else if(pPDG>-3335 && pPDG<-1110) // ==> AntiBaryons
711  {
712  if(oxs) res=aBaCSmanager->GetHMaxT(); // Calculate the max_t value
713  else res=aBaCSmanager->GetExchangeT(Z, N, pPDG);// functionallyRandomized -t in MeV^2
714  }
715  else G4cout<<"*Warning*G4QElastic::CalculateXSt: (i) wrong PDG="<<pPDG<<G4endl;
716  }
717  else G4cout<<"*Warning*G4QElastic::CalculateXSt:*NotInitiatedScattering*"<<G4endl;
718  return res;
719 }