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

Last change on this file since 898 was 819, checked in by garnier, 16 years ago

import all except CVS

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