Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4PiMinusAbsorptionAtRest.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 // File name: G4PiMinusAbsorptionAtRest.hh
27 //
28 // Author: Maria Grazia Pia (pia@genova.infn.it)
29 //
30 // Creation date: 8 May 1998
31 //
32 // Modifications:
33 // MGP 4 Jul 1998 Changed excitation energy calculation
34 // MGP 14 Sep 1998 Fixed excitation energy calculation
35 //
36 // -------------------------------------------------------------------
37 
38 #include "G4ios.hh"
39 
41 
42 #include "G4SystemOfUnits.hh"
43 #include "G4PiMinusStopLi.hh"
44 #include "G4PiMinusStopC.hh"
45 #include "G4PiMinusStopN.hh"
46 #include "G4PiMinusStopO.hh"
47 #include "G4PiMinusStopAl.hh"
48 #include "G4PiMinusStopCu.hh"
49 #include "G4PiMinusStopCo.hh"
50 #include "G4PiMinusStopTa.hh"
51 #include "G4PiMinusStopPb.hh"
54 #include "G4DynamicParticle.hh"
56 #include "Randomize.hh"
57 #include "G4ThreeVector.hh"
58 #include "G4LorentzVector.hh"
60 #include "G4HadronicDeprecate.hh"
61 
62 
63 // Constructor
64 
65 G4PiMinusAbsorptionAtRest::G4PiMinusAbsorptionAtRest(const G4String& processName,
66  G4ProcessType aType) :
67  G4VRestProcess (processName, aType)
68 {
69  G4HadronicDeprecate("G4PiMinusAbsorptionAtRest");
70 
72 
73  _indexDeexcitation = 0;
74 
75  if (verboseLevel>0)
76  { G4cout << GetProcessName() << " is created "<< G4endl; }
78 }
79 
80 
81 // Destructor
82 
84 {
86 }
87 
89 {
91 }
92 
94 {
96 }
97 
99 {
100  const G4DynamicParticle* stoppedHadron = track.GetDynamicParticle();
101 
102  // Check applicability
103  if (! IsApplicable(*(stoppedHadron->GetDefinition())))
104  {
105  G4cerr
106  << "G4PiMinusAbsorptionAtRest: ERROR, particle must be a pion minus!"
107  << G4endl;
108  return NULL;
109  }
110 
111  // Get the current material
112  const G4Material* material = track.GetMaterial();
113 
114  G4double A=-1;
115  G4double Z=-1;
116  G4double random = G4UniformRand();
117  const G4ElementVector* theElementVector = material->GetElementVector();
118  unsigned int i;
119  G4double sum = 0;
120  G4double totalsum=0;
121  for(i=0; i<material->GetNumberOfElements(); ++i)
122  {
123  if((*theElementVector)[i]->GetZ()!=1) totalsum+=material->GetFractionVector()[i];
124  }
125  for (i = 0; i<material->GetNumberOfElements(); ++i)
126  {
127  if((*theElementVector)[i]->GetZ()!=1) sum += material->GetFractionVector()[i];
128  if ( sum/totalsum > random )
129  {
130  A = (*theElementVector)[i]->GetA()*mole/g;
131  Z = (*theElementVector)[i]->GetZ();
132  break;
133  }
134  }
135 
136  // Do the interaction with the nucleon cluster
137 
138  G4PiMinusStopMaterial* algorithm = LoadAlgorithm(static_cast<G4int>(Z));
139  G4PiMinusStopAbsorption stopAbsorption(algorithm,Z,A);
140  stopAbsorption.SetVerboseLevel(verboseLevel);
141 
142  G4DynamicParticleVector* absorptionProducts = stopAbsorption.DoAbsorption();
143 
144  // Deal with the leftover nucleus
145 
146  G4double pionEnergy = stoppedHadron->GetTotalEnergy();
147  G4double excitation = pionEnergy - stopAbsorption.Energy();
148  if (excitation < 0.)
149  {
150  G4Exception("G4PiMinusAbsorptionAtRest::AtRestDoIt()", "HAD_STOP_0000",
151  FatalException, "Excitation energy < 0");
152  }
153  if (verboseLevel>0) { G4cout << " excitation " << excitation << G4endl; }
154 
155  G4StopDeexcitationAlgorithm* nucleusAlgorithm = LoadNucleusAlgorithm();
156  G4StopDeexcitation stopDeexcitation(nucleusAlgorithm);
157 
158  G4double newZ = Z - stopAbsorption.NProtons();
159  G4double newN = A - Z - stopAbsorption.NNeutrons();
160  G4double newA = newZ + newN;
161  G4ReactionProductVector* fragmentationProducts = stopDeexcitation.DoBreakUp(newA,newZ,excitation,stopAbsorption.RecoilMomentum());
162 
163  unsigned int nAbsorptionProducts = 0;
164  if (absorptionProducts != 0)
165  { nAbsorptionProducts = absorptionProducts->size(); }
166 
167  unsigned int nFragmentationProducts = 0;
168  if (fragmentationProducts != 0)
169  { nFragmentationProducts = fragmentationProducts->size(); }
170 
171  if (verboseLevel>0)
172  {
173  G4cout << "nAbsorptionProducts = " << nAbsorptionProducts
174  << " nFragmentationProducts = " << nFragmentationProducts
175  << G4endl;
176  }
177 
178  // Deal with ParticleChange final state
179 
181  aParticleChange.SetNumberOfSecondaries(G4int(nAbsorptionProducts + nFragmentationProducts));
182 
183  for (i = 0; i<nAbsorptionProducts; i++)
184  { aParticleChange.AddSecondary((*absorptionProducts)[i]); }
185 
186 // for (i = 0; i<nFragmentationProducts; i++)
187 // { aParticleChange.AddSecondary(fragmentationProducts->at(i)); }
188  for(i=0; i<nFragmentationProducts; i++)
189  {
190  G4DynamicParticle * aNew =
191  new G4DynamicParticle((*fragmentationProducts)[i]->GetDefinition(),
192  (*fragmentationProducts)[i]->GetMomentum());
193  G4double newTime = aParticleChange.GetGlobalTime((*fragmentationProducts)[i]->GetFormationTime());
194  aParticleChange.AddSecondary(aNew, newTime);
195  delete (*fragmentationProducts)[i];
196  }
197 
198  if (fragmentationProducts != 0) delete fragmentationProducts;
199 
200  if (_indexDeexcitation == 1) aParticleChange.ProposeLocalEnergyDeposit(excitation);
201 
202  // Kill the absorbed pion
204 
205  return &aParticleChange;
206 
207 }
208 
209 G4PiMinusStopMaterial* G4PiMinusAbsorptionAtRest::LoadAlgorithm(int Z)
210 {
211  if (verboseLevel>0)
212  {
213  G4cout << "Load material algorithm " << Z << G4endl;
214  }
215 
216  G4int index = 0;
217  if (Z > 0 && Z < 4) {index = 3;}
218  if (Z > 3 && Z < 7) {index = 6;}
219  if (Z == 7) {index = 7;}
220  if (Z >= 8 && Z<= 11) {index = 8;}
221  if (Z >= 12 && Z<= 18) {index = 13;}
222  if (Z >=19 && Z<= 27) {index = 27;}
223  if (Z >= 28 && Z<= 51) {index = 29;}
224  if (Z >=52 ) {index = 73;}
225 
226  switch (index)
227  {
228  case 3:
229  if (verboseLevel>0)
230  { G4cout << " =================== Load Li algorithm " << G4endl; }
231  return new G4PiMinusStopLi();
232  case 6:
233  if (verboseLevel>0)
234  { G4cout << " =================== Load C algorithm " << G4endl; }
235  return new G4PiMinusStopC();
236  case 7:
237  if (verboseLevel>0)
238  { G4cout << " =================== Load N algorithm " << G4endl; }
239  return new G4PiMinusStopN();
240  case 8:
241  if (verboseLevel>0)
242  { G4cout << " =================== Load O algorithm " << G4endl; }
243  return new G4PiMinusStopO();
244  case 13:
245  if (verboseLevel>0)
246  { G4cout << " =================== Load Al algorithm " << G4endl; }
247  return new G4PiMinusStopAl();
248  case 27:
249  if (verboseLevel>0)
250  { G4cout << " =================== Load Cu algorithm " << G4endl; }
251  return new G4PiMinusStopCu();
252  case 29:
253  if (verboseLevel>0)
254  { G4cout << " =================== Load Co algorithm " << G4endl; }
255  return new G4PiMinusStopCo();
256  case 73:
257  if (verboseLevel>0)
258  { G4cout << " =================== Load Ta algorithm " << G4endl; }
259  return new G4PiMinusStopTa();
260  default:
261  if (verboseLevel>0)
262  { G4cout << " =================== Load default material algorithm " << G4endl; }
263  return new G4PiMinusStopC();
264  }
265 }
266 
267 G4StopDeexcitationAlgorithm* G4PiMinusAbsorptionAtRest::LoadNucleusAlgorithm()
268 {
269 
270  switch (_indexDeexcitation)
271  {
272  case 0:
273  if (verboseLevel>0)
274  { G4cout << " =================== Load Theo deexcitation " << G4endl; }
275  return new G4StopTheoDeexcitation();
276  case 1:
277  if (verboseLevel>0)
278  { G4cout << " =================== Load Dummy deexcitation " << G4endl; }
279  return new G4StopDummyDeexcitation();
280  default:
281  if (verboseLevel>0)
282  { G4cout << " =================== Load default deexcitation " << G4endl; }
283  return new G4StopTheoDeexcitation();
284  }
285 }
286 
288 {
289  _indexDeexcitation = index;
290 }