Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4QNGamma.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 // ---------------- G4QNGamma class -----------------
29 // by Mikhail Kossov, December 2003.
30 // G4QNGamma class of the CHIPS Simulation Branch in GEANT4
31 // ---------------------------------------------------------------
32 // **************************************************************************
33 // This Header is a part of the CHIPS Physics Package (author: M. Kosov)
34 // **************************************************************************
35 // Short description: This is a universal class for the incoherent (inelastic)
36 // nuclear (n,gamma) interactions (neutron capture) in the CHIPS model.
37 // @@ At present the gamma cascade is not simulated (one final photon)
38 // ---------------------------------------------------------------------------
39 //#define debug
40 //#define pdebug
41 
42 #include "G4QNGamma.hh"
43 #include "G4HadronicDeprecate.hh"
44 
45 
46 // Initialization of static Material/Element/Isotope vectors
47 std::vector<G4int> G4QNGamma::ElementZ; // Z of the element(i) in the Last Calc
48 std::vector<G4double> G4QNGamma::ElProbInMat; // ProbabilitySum ofElements inMaterial
49 std::vector<std::vector<G4int>*> G4QNGamma::ElIsoN;// # of isotope(j), # of Element(i)
50 std::vector<std::vector<G4double>*>G4QNGamma::IsoProbInEl;//SumProbabIsotopes in Element I
51 
52 // Constructor
53 G4QNGamma::G4QNGamma(const G4String& processName)
54  : G4VDiscreteProcess(processName, fHadronic)
55 {
56  G4HadronicDeprecate("G4QNGamma");
57 
58  EnMomConservation = G4LorentzVector(0.,0.,0.,0.);
59  nOfNeutrons = 0;
60 #ifdef debug
61  G4cout<<"G4QNGamma::Constructor is called"<<G4endl;
62 #endif
63  if (verboseLevel>0) G4cout << GetProcessName() << " process is created "<< G4endl;
64 }
65 
66 // Destructor (standard procedure @@ to be moved to G4VQProcess)
68 {
69  // The following is just a copy of what is done in PostStepDoIt every interaction !
70  // The correction is if(IPIE), so just for(...;ip<IPIE;...) does not work ! @@
71  G4int IPIE=IsoProbInEl.size(); // How many old elements?
72  if(IPIE) for(G4int ip=0; ip<IPIE; ++ip) // Clean up the SumProb's of Isotopes (SPI)
73  {
74  std::vector<G4double>* SPI=IsoProbInEl[ip]; // Pointer to the SPI vector
75  SPI->clear();
76  delete SPI;
77  std::vector<G4int>* IsN=ElIsoN[ip]; // Pointer to the N vector
78  IsN->clear();
79  delete IsN;
80  }
81  ElProbInMat.clear(); // Clean up the SumProb's of Elements (SPE)
82  ElementZ.clear(); // Clear the body vector for Z of Elements
83  IsoProbInEl.clear(); // Clear the body vector for SPI
84  ElIsoN.clear(); // Clear the body vector for N of Isotopes
85 }
86 
87 
89 {
90  return EnMomConservation;
91 }
92 
93 G4int G4QNGamma::GetNumberOfNeutronsInTarget() // @@ move to G4VQProcess
94 {
95  return nOfNeutrons;
96 }
97 
98 // output of the function must be in units of length! L=1/sig_V,sig_V=SUM(n(j,i)*sig(j,i)),
99 // where n(i,j) is a number of nuclei of the isotop j of the element i in V=1(lengtUnit^3)
100 // ********** All CHIPS cross sections are calculated in the surface units ************
101 // @@ Can demand 3 internal functions when G4VQProcess is used @@ future plans
103 {
104 #ifdef debug
105  G4cout<<"G4QNGamma::GetMeanFreePath: Called Fc="<<*Fc<<G4endl;
106 #endif
107  *Fc = NotForced;
108 #ifdef debug
109  G4cout<<"G4QNGamma::GetMeanFreePath: Before GetDynPart"<<G4endl;
110 #endif
111  const G4DynamicParticle* incidentParticle = aTrack.GetDynamicParticle();
112 #ifdef debug
113  G4cout<<"G4QNGamma::GetMeanFreePath: Before GetDef"<<G4endl;
114 #endif
115  G4ParticleDefinition* incidentParticleDefinition=incidentParticle->GetDefinition();
116  G4double Momentum = incidentParticle->GetTotalMomentum(); // 3-momentum of the Particle
117  if( !IsApplicable(*incidentParticleDefinition)) // @@ Unique for all QProcesses
118  {
119  G4cout<<"-W-G4QNGamma::GetMeanFreePath called for not implemented particle"<<G4endl;
120  return DBL_MAX;
121  }
122 #ifdef debug
123  G4cout<<"G4QNGamma::GetMeanFreePath: BeforeGetMaterial P="<<Momentum<<G4endl;
124 #endif
125  // @@ Can be additional condition internal function of G4VQProcess
126  if(Momentum > 500.) return DBL_MAX; // @@ Temporary cut (QInternal=MeV -> IU!)
127  // @@ This is a standard procedure, which can be moved to G4VQProcess (above is a funct)
128  const G4Material* material = aTrack.GetMaterial(); // Get the current material
129  const G4double* NOfNucPerVolume = material->GetVecNbOfAtomsPerVolume();
130  const G4ElementVector* theElementVector = material->GetElementVector();
131  G4int nE=material->GetNumberOfElements();
132 #ifdef debug
133  G4cout<<"G4QNGamma::GetMeanFreePath:"<<nE<<" Elem's in theMaterial"<<G4endl;
134 #endif
135  // @@ Can be internal function called by GetMeanFreePath (above Isotope LOOP)
136  G4VQCrossSection* CSmanager = 0; // @@ Reference modified in the function
137  G4QNeutronCaptureRatio* capMan = 0; // @@ Reference modified in the function
138  G4int pPDG =0; // @@ Reference modified in the function
139  G4double sigma =0.; // CS mean over isotopes @@ Reference modified in the function
140  if(incidentParticleDefinition == G4Neutron::Neutron())
141  {
143  capMan=G4QNeutronCaptureRatio::GetPointer(); // @@ can be CSmanager2
144 #ifdef debug
145  G4cout<<"G4QNGamma::GetMeanFreePath: CSmanager is defined for neutrons"<<G4endl;
146 #endif
147  pPDG=2112;
148  }
149  else
150  {
151  G4cout<<"-Warning-G4QNGamma::GetMeanFreePath:Particle "
152  <<incidentParticleDefinition->GetPDGEncoding()<<" isn't a neutron"<<G4endl;
153  return DBL_MAX; // can be returned in sigma
154  }
155  // @@ End of possible internal function
156  G4QIsotope* Isotopes = G4QIsotope::Get(); // Pointer to the G4QIsotopes singleton
157  G4int IPIE=IsoProbInEl.size(); // How many old elements?
158  if(IPIE) for(G4int ip=0; ip<IPIE; ++ip) // Clean up the SumProb's of Isotopes (SPI)
159  {
160  std::vector<G4double>* SPI=IsoProbInEl[ip]; // Pointer to the SPI vector
161  SPI->clear();
162  delete SPI;
163  std::vector<G4int>* IsN=ElIsoN[ip]; // Pointer to the N vector
164  IsN->clear();
165  delete IsN;
166  }
167  ElProbInMat.clear(); // Clean up the SumProb's of Elements (SPE)
168  ElementZ.clear(); // Clear the body vector for Z of Elements
169  IsoProbInEl.clear(); // Clear the body vector for SPI
170  ElIsoN.clear(); // Clear the body vector for N of Isotopes
171  for(G4int i=0; i<nE; ++i)
172  {
173  G4Element* pElement=(*theElementVector)[i]; // Pointer to the current element
174  G4int Z = static_cast<G4int>(pElement->GetZ()); // Z of the Element
175  ElementZ.push_back(Z); // Remember Z of the Element
176  G4int isoSize=0; // The default for the isoVectorLength is 0
177  G4int indEl=0; // Index of non-trivial element or 0(default)
178  G4IsotopeVector* isoVector=pElement->GetIsotopeVector(); // Get the predefined IsoVect
179  if(isoVector) isoSize=isoVector->size();// Get size of the existing isotopeVector
180 #ifdef debug
181  G4cout<<"G4QNGamma::GetMeanFreePath: isovectorLength="<<isoSize<<G4endl; // Result
182 #endif
183  if(isoSize) // The Element has non-trivial abundance set
184  {
185  indEl=pElement->GetIndex()+1; // Index of the non-trivial element
186  if(!Isotopes->IsDefined(Z,indEl)) // This index is not defined for this Z: define
187  {
188  std::vector<std::pair<G4int,G4double>*>* newAbund =
189  new std::vector<std::pair<G4int,G4double>*>;
190  G4double* abuVector=pElement->GetRelativeAbundanceVector();
191  for(G4int j=0; j<isoSize; j++) // Calculation of abundance vector for isotopes
192  {
193  G4int N=pElement->GetIsotope(j)->GetN()-Z; // N means A=N+Z !
194  if(pElement->GetIsotope(j)->GetZ()!=Z)G4cerr<<"G4QNGamma::GetMeanFreePath"
195  <<": Z="<<pElement->GetIsotope(j)->GetZ()<<"#"<<Z<<G4endl;
196  G4double abund=abuVector[j];
197  std::pair<G4int,G4double>* pr= new std::pair<G4int,G4double>(N,abund);
198 #ifdef debug
199  G4cout<<"G4QNGamma::GetMeanFreePath: p#="<<j<<",N="<<N<<",ab="<<abund<<G4endl;
200 #endif
201  newAbund->push_back(pr);
202  }
203 #ifdef debug
204  G4cout<<"G4QNGamma::GetMeanFreePath: pairVectLength="<<newAbund->size()<<G4endl;
205 #endif
206  indEl=G4QIsotope::Get()->InitElement(Z,indEl,newAbund); // definition of the newInd
207  for(G4int k=0; k<isoSize; k++) delete (*newAbund)[k]; // Cleaning temporary
208  delete newAbund; // Was "new" in the beginning of the name space
209  }
210  }
211  std::vector<std::pair<G4int,G4double>*>* cs= Isotopes->GetCSVector(Z,indEl);//CSPointer
212  std::vector<G4double>* SPI = new std::vector<G4double>; // Pointer to the SPI vector
213  IsoProbInEl.push_back(SPI);
214  std::vector<G4int>* IsN = new std::vector<G4int>; // Pointer to the N vector
215  ElIsoN.push_back(IsN);
216  G4int nIs=cs->size(); // A#Of Isotopes in the Element
217  G4double susi=0.; // sum of CS over isotopes
218 #ifdef debug
219  G4cout<<"G4QNGamma::GetMeanFreePath: Before Loop nIs="<<nIs<<G4endl;
220 #endif
221  if(nIs) for(G4int j=0; j<nIs; j++) // Calculate CS for eachIsotope of El
222  {
223  std::pair<G4int,G4double>* curIs=(*cs)[j]; // A pointer, which is used twice
224  G4int N=curIs->first; // #of Neuterons in the isotope j of El i
225  IsN->push_back(N); // Remember Min N for the Element
226 #ifdef debug
227  G4cout<<"G4QNGam::GetMeanFrP: Before CS, P="<<Momentum<<",Z="<<Z<<",N="<<N<<G4endl;
228 #endif
229  // @@ Can be a function, depanding on CSm1, CSm2, Momentum, Z, N, pPDG
230  G4double CSI=CSmanager->GetCrossSection(true, Momentum, Z, N, pPDG) *
231  capMan->GetRatio(Momentum, Z, N); // CS(j,i) for the isotope
232 #ifdef debug
233  G4cout<<"GQC::GMF:X="<<CSI<<",M="<<Momentum<<",Z="<<Z<<",N="<<N<<",P="<<pPDG<<G4endl;
234 #endif
235  curIs->second = CSI;
236  susi+=CSI; // Make a sum per isotopes
237  SPI->push_back(susi); // Remember summed cross-section
238  } // End of temporary initialization of the cross sections in the G4QIsotope singeltone
239  sigma+=Isotopes->GetMeanCrossSection(Z,indEl)*NOfNucPerVolume[i];//SUM(MeanCS*NOfNperV)
240  ElProbInMat.push_back(sigma);
241  } // End of LOOP over Elements
242 #ifdef debug
243  G4cout<<"G4QNGam::GetMeanFrPa: Sigma="<<sigma<<G4endl;
244 #endif
245  if(sigma > 0.) return 1./sigma; // Mean path [distance]
246  return DBL_MAX;
247 }
248 
249 // Original in any G4QProcess inheriting G4Process
251 {
252  if ( particle == *( G4Neutron::Neutron() ) ) return true;
253 #ifdef debug
254  G4cout<<"***G4QNGamma::IsApplicable: PDG="<<particle.GetPDGEncoding()<<G4endl;
255 #endif
256  return false;
257 }
258 
260 {
261 #ifdef debug
262  static const G4double mNeut= G4QPDGCode(2112).GetMass();
263 #endif
264  static const G4LorentzVector vacuum4M(0.,0.,0.,0.);
265  //-------------------------------------------------------------------------------------
266  const G4DynamicParticle* projHadron = track.GetDynamicParticle();
267  const G4ParticleDefinition* particle=projHadron->GetDefinition();
268 #ifdef debug
269  G4cout<<"G4QNGamma::PostStepDoIt: Before the GetMeanFreePath is called"<<G4endl;
270 #endif
272  GetMeanFreePath(track, 1., &cond);
273 #ifdef debug
274  G4cout<<"G4QNGamma::PostStepDoIt: After the GetMeanFreePath is called"<<G4endl;
275 #endif
276  G4LorentzVector proj4M=projHadron->Get4Momentum(); // 4-momentum of the projectile (IU?)
277  G4double momentum = projHadron->GetTotalMomentum(); // 3-momentum of the Particle
278  G4double Momentum=proj4M.rho();
279  if(std::fabs(Momentum-momentum)>.001)
280  G4cerr<<"*G4QNGamma::PostStepDoIt: P="<<Momentum<<"#"<<momentum<<G4endl;
281 #ifdef debug
282  G4double mp=proj4M.m(); // @@ must be just the neutron mass
283  if(std::fabs(mp-mNeut)>.001)G4cerr<<"*G4QNGamma::PostStDoIt: M="<<mp<<"#"<<mNeut<<G4endl;
284  G4cout<<"->G4QNGam::PostStDoIt:*called*,4M="<<proj4M<<",P="<<Momentum<<",m="<<mp<<G4endl;
285 #endif
286  // The same cut function can be used as in MeanFreePath (500)
287  if (!IsApplicable(*particle) || Momentum > 500.) // Check applicability (@@ IU?)
288  {
289  G4cerr<<"G4QNGamma::PostStepDoIt: Only neutrons with P="<<Momentum<<" < 500"<<G4endl;
290  return 0;
291  }
292  const G4Material* material = track.GetMaterial(); // Get the current material
293  G4int Z=0;
294  const G4ElementVector* theElementVector = material->GetElementVector();
295  G4int nE=material->GetNumberOfElements();
296 #ifdef debug
297  G4cout<<"G4QNGamma::PostStepDoIt: "<<nE<<" elements in the material."<<G4endl;
298 #endif
299  G4int EPIM=ElProbInMat.size();
300 #ifdef debug
301  G4cout<<"G4QNGam::PostStDoIt: m="<<EPIM<<",n="<<nE<<",T="<<ElProbInMat[EPIM-1]<<G4endl;
302 #endif
303  G4int i=0;
304  if(EPIM>1)
305  {
306  G4double rnd = ElProbInMat[EPIM-1]*G4UniformRand();
307  for(i=0; i<nE; ++i)
308  {
309 #ifdef debug
310  G4cout<<"G4QNGamma::PostStepDoIt:E["<<i<<"]="<<ElProbInMat[i]<<",r="<<rnd<<G4endl;
311 #endif
312  if (rnd<ElProbInMat[i]) break;
313  }
314  if(i>=nE) i=nE-1; // Top limit for the Element
315  }
316  G4Element* pElement=(*theElementVector)[i];
317  Z=static_cast<G4int>(pElement->GetZ());
318 #ifdef debug
319  G4cout<<"G4QNGamma::PostStepDoIt: i="<<i<<", Z(element)="<<Z<<G4endl;
320 #endif
321  if(Z <= 0)
322  {
323  G4cerr<<"---Warning---G4QNGamma::PostStepDoIt: Element with Z="<<Z<<G4endl;
324  if(Z<0) return 0;
325  }
326  std::vector<G4double>* SPI = IsoProbInEl[i];// Vector of summedProbabilities for isotopes
327  std::vector<G4int>* IsN = ElIsoN[i]; // Vector of "#of neutrons" in the isotope El[i]
328  G4int nofIsot=SPI->size(); // #of isotopes in the element i
329 #ifdef debug
330  G4cout<<"G4QNGam::PosStDoIt:n="<<nofIsot<<",T="<<(*SPI)[nofIsot-1]<<G4endl;
331 #endif
332  G4int j=0;
333  if(nofIsot>1)
334  {
335  G4double rndI=(*SPI)[nofIsot-1]*G4UniformRand(); // Randomize the isotop of the Element
336  for(j=0; j<nofIsot; ++j)
337  {
338 #ifdef debug
339  G4cout<<"G4QNGamma::PostStepDoIt: SP["<<j<<"]="<<(*SPI)[j]<<", r="<<rndI<<G4endl;
340 #endif
341  if(rndI < (*SPI)[j]) break;
342  }
343  if(j>=nofIsot) j=nofIsot-1; // Top limit for the isotope
344  }
345  G4int N =(*IsN)[j]; // Randomized number of neutrons
346 #ifdef debug
347  G4cout<<"G4QNGamma::PostStepDoIt: Z="<<Z<<", j="<<i<<", N(isotope)="<<N<<G4endl;
348 #endif
349  G4double kinEnergy= projHadron->GetKineticEnergy();
351  //if() //DoNothing Action insead of the reaction
352  //{
353  // aParticleChange.ProposeEnergy(kinEnergy);
354  // aParticleChange.ProposeLocalEnergyDeposit(0.);
355  // aParticleChange.ProposeMomentumDirection(dir);
356  // aParticleChange.ProposeTrackStatus(fAlive);
357  // return G4VDiscreteProcess::PostStepDoIt(track,step);
358  //}
359  if(N<0)
360  {
361  G4cerr<<"-Warning-G4QNGamma::PostStepDoIt: Isotope with Z="<<Z<<", 0>N="<<N<<G4endl;
362  return 0;
363  }
364  nOfNeutrons=N; // Remember it for the energy-momentum check
365 #ifdef debug
366  G4cout<<"G4QNGamma::PostStepDoIt: N="<<N<<" for element with Z="<<Z<<G4endl;
367 #endif
369  G4double weight = track.GetWeight();
370 #ifdef debug
371  G4cout<<"G4QNGamma::PostStepDoIt: weight="<<weight<<G4endl;
372 #endif
373  G4double localtime = track.GetGlobalTime();
374 #ifdef debug
375  G4cout<<"G4QNGamma::PostStepDoIt: localtime="<<localtime<<G4endl;
376 #endif
378  G4TouchableHandle trTouchable = track.GetTouchableHandle();
379 #ifdef debug
380  G4cout<<"G4QNGamma::PostStepDoIt: position="<<position<<G4endl;
381 #endif
382  G4int targPDG = 90000000 + Z*1000 + N; // PDG Code of the target nucleus
383  G4QPDGCode targQPDG(targPDG);
384  G4double tM = targQPDG.GetMass(); // Target mass
385 #ifdef debug
386  G4cout<<"G4QNGamma::PostStepDoIt: n + targPDG="<<targPDG<<G4endl;
387 #endif
388  // @@ All above is universal for all processes except for the additional condition (500)
389  G4LorentzVector tot4M=G4LorentzVector(0.,0.,0.,tM)+proj4M;
390  G4double totM2=tot4M.m2();
391  G4int tZ=Z;
392  G4int tN=N+1;
393  G4int resPDG = targPDG + 1; // Final ++N nucleus PDG
394  G4double rM=G4QPDGCode(resPDG).GetMass(); // Mass of the final nucleus
395  G4LorentzVector r4M=G4LorentzVector(0.,0.,0.,rM); // 4mom of the final nucleus
396  G4LorentzVector g4M=G4LorentzVector(0.,0.,0.,0.); // 4mom of the gamma
397 #ifdef debug
398  G4cout<<"G4QNGamma::PostStepDoIt: tM="<<tM << ", rM="<<rM << ", Q="<<tM+mNeut-rM<<G4endl;
399 #endif
400  if(!G4QHadron(tot4M).DecayIn2(r4M, g4M)) // The compoun decay din't succeed
401  {
402  //G4cerr<<"G4QNGamma::PostStDoIt: tM="<<std::sqrt(totM2)<<" < rM="<<rM<<G4endl;
403  //G4Exception("G4QNGamma::PostStepDoIt()", "HAD_CHPS_0001",
404  // FatalException, "Hadronize quasmon: Can't decay TotNuc->ResNuc+gam");
405  G4cerr<<"-Warning-G4QNGamma::PostStDoIt: tM="<<std::sqrt(totM2)<<" < rM="<<rM<<G4endl;
406  aParticleChange.ProposeEnergy(kinEnergy);
410  return G4VDiscreteProcess::PostStepDoIt(track,step);
411  }
412 #ifdef debug
413  G4cout<<"G4QNGam::PStDoIt: RA="<<r4M.rho()<<r4M<<", Gamma="<<g4M.rho()<<g4M<<G4endl;
414 #endif
415  EnMomConservation = tot4M - r4M - g4M; // EM conservation check 4mom
416  aParticleChange.ProposeEnergy(0.); // A standard procedure of killing proj.
417  aParticleChange.ProposeTrackStatus(fStopAndKill);// projectile neutron is killed
418  aParticleChange.SetNumberOfSecondaries(2); // Fix a#of secondaries
419  // Fill the gamma
420  G4ParticleDefinition* theDefinition = G4Gamma::Gamma();
421  G4DynamicParticle* theGam = new G4DynamicParticle(theDefinition, g4M);
422  G4Track* capGamma = new G4Track(theGam, localtime, position );
423  capGamma->SetWeight(weight);
424  capGamma->SetTouchableHandle(trTouchable);
425  aParticleChange.AddSecondary(capGamma);
426  // ----------------------------------------------------
427  // Fill the final nucleus
428  G4int tA=tZ+tN;
429  if (resPDG==90000001) theDefinition = G4Neutron::Neutron();
430  else if(resPDG==90001000) theDefinition = G4Proton::Proton();
431  else theDefinition = G4ParticleTable::GetParticleTable()->FindIon(tZ, tA, 0, tZ);
432  G4DynamicParticle* theReN = new G4DynamicParticle(theDefinition, r4M);
433  G4Track* scatReN = new G4Track(theReN, localtime, position );
434  scatReN->SetWeight(weight);
435  scatReN->SetTouchableHandle(trTouchable);
436  aParticleChange.AddSecondary(scatReN);
437 
438  return G4VDiscreteProcess::PostStepDoIt(track, step);
439 }