Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4PiMinusStopCo.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: G4PiMinusStopCo
27//
28// Author: Maria Grazia Pia (pia@genova.infn.it)
29//
30// Creation date: 8 May 1998
31//
32// -------------------------------------------------------------------
33
34#include "G4ios.hh"
35
36#include "G4PiMinusStopCo.hh"
37
38#include <vector>
39
40#include "globals.hh"
41#include "Randomize.hh"
42#include "G4Proton.hh"
43#include "G4Neutron.hh"
44#include "G4ParticleTypes.hh"
47#include "G4LorentzVector.hh"
50
51// np/pp production ratio
52// Experimental values:
53// R(np/pp) (E. Gadioli et al., Phys Rev C36 (1987) 741
54G4double G4PiMinusStopCo::npRatio = 2.5;
55
56
57
58// Average numbers of final nucleons detected, for N-pair absorption
59// (Hartmann et al., Nucl. Phys. A300 (1978) 345
60G4double G4PiMinusStopCo::nFinalNucleons = 1.38;
61
62// Kinetic energy (MeV) distributions measured for coincident nucleon
63// emission
64// P. Heusi et al., Nucl. Phys. A407(1983) 429
65
66G4int G4PiMinusStopCo::eKinEntries = 11;
67
68G4double G4PiMinusStopCo::eKinData[11] = {0.085, 0.09, 0.09, 0.15,
69 0.1, 0.09, 0.08,
70 0.04, 0.03, 0.02, 0.01};
71
72G4double G4PiMinusStopCo::eKin[12] = { 15., 17.5, 25., 33.,
73 42., 52., 62.,
74 75., 85., 95., 105. };
75
76
77// Opening angle distributions measured for coincident nucleon emission
78// (P.Heusi et al., Nucl. Phys. A407 (1983) 429
79
80G4int G4PiMinusStopCo::angleEntries = 7;
81
82G4double G4PiMinusStopCo::angleData[7] =
83{6., 8., 9., 10., 25., 40., 45. };
84
85G4double G4PiMinusStopCo::angle[8] = { 0.5, 0.7, 0.87, 1.4, 2.1,
86 2.44, 2.8, 3.1415927 };
87
88
89
90// Constructor
91
93
94{
95 // Cluster size: nucleon pair, alpha, triton etc.
96 // First implementation: interaction with nucleon pair only
97 _clusterSize = 2;
98
99 // R ratio
100 theR = 1. / (1. + npRatio);
101
102 _definitions = new std::vector<G4ParticleDefinition*>();
103 _momenta = new std::vector<G4LorentzVector*>();
104
105 std::vector<double> eKinVector;
106 std::vector<double> eKinDataVector;
107 int i;
108 for (i=0; i<eKinEntries; i++)
109 {
110 eKinVector.push_back(eKin[i]);
111 eKinDataVector.push_back(eKinData[i]);
112 }
113 eKinVector.push_back(eKin[eKinEntries]);
114 _distributionE = new G4DistributionGenerator(eKinVector,eKinDataVector);
115
116 std::vector<double> angleVector;
117 std::vector<double> angleDataVector;
118 for (i=0; i<angleEntries; i++)
119 {
120 angleVector.push_back(angle[i]);
121 angleDataVector.push_back(angleData[i]);
122 }
123 angleVector.push_back(angle[angleEntries]);
124 _distributionAngle = new G4DistributionGenerator(angleVector,angleDataVector);
125}
126
127
128// Destructor
129
131{}
132
134{
135 return nFinalNucleons;
136}
137
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
virtual G4double FinalNucleons()
virtual ~G4PiMinusStopCo()
std::vector< G4ParticleDefinition * > * _definitions
G4DistributionGenerator * _distributionE
std::vector< G4LorentzVector * > * _momenta
G4DistributionGenerator * _distributionAngle