source: trunk/source/processes/hadronic/stopping/src/G4PiMinusStopAbsorption.cc @ 1307

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

update processes

File size: 6.7 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//      File name:     G4PiMinusStopAbsorption
27//
28//      Author:        Maria Grazia Pia (pia@genova.infn.it)
29//
30//      Creation date: 8 May 1998
31//
32// -------------------------------------------------------------------
33
34
35#include "G4PiMinusStopAbsorption.hh"
36#include <vector>
37
38#include "globals.hh"
39#include "Randomize.hh"
40#include "G4NucleiProperties.hh"
41#include "G4ParticleTypes.hh"
42#include "G4Nucleus.hh"
43#include "G4ReactionKinematics.hh"
44#include "G4DynamicParticleVector.hh"
45#include "G4ParticleDefinition.hh"
46#include "G4Proton.hh"
47#include "G4Neutron.hh"
48#include "G4ThreeVector.hh"
49
50// Constructor
51
52G4PiMinusStopAbsorption::G4PiMinusStopAbsorption(G4PiMinusStopMaterial* materialAlgo, 
53                                                 const G4double Z, const G4double A)
54 
55{
56  _materialAlgo = materialAlgo;
57  _nucleusZ = Z;
58  _nucleusA = A;
59  _level = 0;
60  _absorptionProducts = new G4DynamicParticleVector();
61}
62
63// Destructor
64
65G4PiMinusStopAbsorption::~G4PiMinusStopAbsorption()
66{
67  // Memory management of materialAlgo needs better thought (MGP)
68  delete _materialAlgo;
69  // Who owns it? Memory management is not clear... (MGP)
70  //  _absorptionProducts->clearAndDestroy();
71  delete _absorptionProducts;
72}
73
74G4DynamicParticleVector* G4PiMinusStopAbsorption::DoAbsorption()
75{
76  std::vector<G4ParticleDefinition*>* defNucleons = _materialAlgo->DefinitionVector();
77
78  G4double newA = _nucleusA;
79  G4double newZ = _nucleusZ;
80
81  if (defNucleons != 0)
82    {
83      for (unsigned int i=0; i<defNucleons->size(); i++)
84        {
85          if ( (*defNucleons)[i] == G4Proton::Proton())
86            {
87              newA = newA - 1;
88              newZ = newZ - 1;
89            }
90          if ((*defNucleons)[i] == G4Neutron::Neutron()) 
91            { newA = newA - 1; }
92        }
93    }
94
95  G4double binding = G4NucleiProperties::GetBindingEnergy(static_cast<G4int>(_nucleusA) ,static_cast<G4int>(_nucleusZ)) / _nucleusA;
96  G4double mass = G4NucleiProperties::GetNuclearMass(static_cast<G4int>(newA),static_cast<G4int>(newZ));
97
98
99  std::vector<G4LorentzVector*>* p4Nucleons = _materialAlgo->P4Vector(binding,mass);
100
101  if (defNucleons != 0 && p4Nucleons != 0)
102    {
103      unsigned int nNucleons = p4Nucleons->size();
104     
105      G4double seen = _materialAlgo->FinalNucleons() / 2.;
106      G4int maxN = nNucleons;
107      if (defNucleons->size() < nNucleons) { maxN = defNucleons->size(); }
108
109      for (G4int i=0; i<maxN; i++)
110        {
111          G4DynamicParticle* product;
112          if ((*defNucleons)[i] == G4Proton::Proton()) 
113            { product = new G4DynamicParticle(G4Proton::Proton(),*((*p4Nucleons)[i])); }
114          else
115            { product = new G4DynamicParticle(G4Neutron::Neutron(),*((*p4Nucleons)[i])); }
116          G4double ranflat = G4UniformRand();
117         
118          if (ranflat < seen) 
119            { _absorptionProducts->push_back(product); }
120          else
121            { delete product; }
122        }
123    }
124
125  return _absorptionProducts;
126
127}
128
129G4ThreeVector G4PiMinusStopAbsorption::RecoilMomentum()
130{
131  G4ThreeVector pProducts(0.,0.,0.);
132 
133  for (unsigned int i = 0; i< _absorptionProducts->size(); i++)
134    {
135      pProducts = pProducts + (*_absorptionProducts)[i]->GetMomentum();
136    }
137  return pProducts;
138}
139
140
141G4int G4PiMinusStopAbsorption::NProtons()
142{
143  G4int n = 0;
144  G4int entries = _absorptionProducts->size();
145  for (int i = 0; i<entries; i++)
146    {
147      if ((*_absorptionProducts)[i]->GetDefinition() == G4Proton::Proton())
148        { n = n + 1; }
149    }
150  return n;
151}
152
153
154G4int G4PiMinusStopAbsorption::NNeutrons()
155{
156  G4int n = 0;
157  G4int entries = _absorptionProducts->size();
158  for (int i = 0; i<entries; i++)
159    {
160      if ((*_absorptionProducts)[i]->GetDefinition() == G4Neutron::Neutron())
161        { n = n + 1; }
162    }
163  return n;
164}
165
166
167G4double G4PiMinusStopAbsorption::Energy()
168{
169  G4double energy = 0.;
170  G4double productEnergy = 0.;
171  G4ThreeVector pProducts(0.,0.,0.);
172  G4int nN = 0;
173  G4int nP = 0;
174
175
176  G4int nAbsorptionProducts = _absorptionProducts->size();
177
178  for (int i = 0; i<nAbsorptionProducts; i++)
179    {
180      productEnergy += (*_absorptionProducts)[i]->GetKineticEnergy();
181      pProducts = pProducts + (*_absorptionProducts)[i]->GetMomentum();
182      if ((*_absorptionProducts)[i]->GetDefinition() == G4Neutron::Neutron()) nN++;
183      if ((*_absorptionProducts)[i]->GetDefinition() == G4Proton::Proton()) nP++;
184    }
185
186  G4double productBinding = (G4NucleiProperties::GetBindingEnergy(static_cast<G4int>(_nucleusA),static_cast<G4int>(_nucleusZ)) / _nucleusA) * nAbsorptionProducts;
187  G4double mass = G4NucleiProperties::GetNuclearMass(_nucleusA - (nP + nN),_nucleusZ - nP);
188  G4double pNucleus = pProducts.mag();
189  G4double eNucleus = std::sqrt(pNucleus*pNucleus + mass*mass);
190  G4double tNucleus = eNucleus - mass;
191  G4double temp = 
192    G4NucleiProperties::GetBindingEnergy(static_cast<G4int>(_nucleusA - (nP + nN)),static_cast<G4int>(_nucleusZ - nP)) - 
193    G4NucleiProperties::GetBindingEnergy(static_cast<G4int>(_nucleusA),static_cast<G4int>(_nucleusZ));
194  energy = productEnergy + productBinding + tNucleus;
195 
196  if (_level > 0)
197    {
198      std::cout << "E products " <<  productEnergy 
199           << " Binding " << productBinding << " " << temp << " "
200           << " Tnucleus " << tNucleus
201           << " energy = " << energy << G4endl;
202    }
203
204  return energy;
205}
206
207void G4PiMinusStopAbsorption::SetVerboseLevel(G4int level)
208{
209  _level = level;
210  return;
211}
Note: See TracBrowser for help on using the repository browser.