Geant4  10.02.p01
G4TwoBodyAngularDist.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: G4TwoBodyAngularDist.cc 71652 2013-06-19 17:20:45Z mkelsey $
27 // Author: Michael Kelsey (SLAC)
28 // Date: 21 February 2013
29 //
30 // Description: Singleton class to evaluate two-body angular distribution
31 // functions based on intial/final state codes.
32 //
33 // 20130307 M. Kelsey -- Add verbosity interface
34 // 20130422 M. Kelsey -- Add three-body distributions, for temporary use
35 // 20130619 Change singleton instance to be thread-local, to avoid collisions.
36 // 20141121 Use G4AutoDelete to avoid end-of-thread memory leaks
37 
38 #include "G4TwoBodyAngularDist.hh"
39 #include "G4AutoDelete.hh"
40 #include "G4GamP2NPipAngDst.hh"
41 #include "G4GamP2PPi0AngDst.hh"
42 #include "G4GammaNuclAngDst.hh"
43 #include "G4PP2PPAngDst.hh"
44 #include "G4NP2NPAngDst.hh"
45 #include "G4Pi0P2Pi0PAngDst.hh"
46 #include "G4PimP2Pi0NAngDst.hh"
47 #include "G4PimP2PimPAngDst.hh"
48 #include "G4PipP2PipPAngDst.hh"
49 #include "G4HadNElastic1AngDst.hh"
50 #include "G4HadNElastic2AngDst.hh"
51 #include "G4InuclParticleNames.hh"
52 #include "G4NuclNuclAngDst.hh"
53 #include "G4HadNucl3BodyAngDst.hh"
54 #include "G4NuclNucl3BodyAngDst.hh"
55 #include "G4PiNInelasticAngDst.hh"
56 #include "G4InuclParticleNames.hh"
57 using namespace G4InuclParticleNames;
58 
59 // Singleton is created at first invocation
60 
62 
64  if (!theInstance) {
65  theInstance = new G4TwoBodyAngularDist;
66  G4AutoDelete::Register(theInstance);
67  }
68 
69  return theInstance;
70 }
71 
72 // Constructor and destructor
73 
75  : gp_npip(new G4GamP2NPipAngDst), gp_ppi0(new G4GamP2PPi0AngDst),
76  ppAngDst(new G4PP2PPAngDst), npAngDst(new G4NP2NPAngDst),
77  nnAngDst(new G4NuclNuclAngDst), pi0pAngDst(new G4Pi0P2Pi0PAngDst),
78  pipCXAngDst(new G4PimP2Pi0NAngDst), pimpAngDst(new G4PimP2PimPAngDst),
79  pippAngDst(new G4PipP2PipPAngDst), qxAngDst(new G4PiNInelasticAngDst),
80  hn1AngDst(new G4HadNElastic1AngDst), hn2AngDst(new G4HadNElastic2AngDst),
81  gnAngDst(new G4GammaNuclAngDst), hn3BodyDst(new G4HadNucl3BodyAngDst),
82  nn3BodyDst(new G4NuclNucl3BodyAngDst)
83 {;}
84 
86  delete gp_npip;
87  delete gp_ppi0;
88  delete ppAngDst;
89  delete nnAngDst;
90  delete pi0pAngDst;
91  delete pipCXAngDst;
92  delete pimpAngDst;
93  delete pippAngDst;
94  delete qxAngDst;
95  delete hn1AngDst;
96  delete hn2AngDst;
97  delete gnAngDst;
98  delete npAngDst;
99  delete hn3BodyDst;
100  delete nn3BodyDst;
101 }
102 
103 
104 // Set verbosity for all generators (const-cast required)
105 
107  const_cast<G4TwoBodyAngularDist*>(GetInstance())->passVerbose(verbose);
108 }
109 
111  if (gp_npip) gp_npip->setVerboseLevel(verbose);
112  if (gp_ppi0) gp_ppi0->setVerboseLevel(verbose);
113  if (ppAngDst) ppAngDst->setVerboseLevel(verbose);
114  if (nnAngDst) nnAngDst->setVerboseLevel(verbose);
115  if (pi0pAngDst) pi0pAngDst->setVerboseLevel(verbose);
116  if (pipCXAngDst) pipCXAngDst->setVerboseLevel(verbose);
117  if (pimpAngDst) pimpAngDst->setVerboseLevel(verbose);
118  if (pippAngDst) pippAngDst->setVerboseLevel(verbose);
119  if (qxAngDst) qxAngDst->setVerboseLevel(verbose);
120  if (hn1AngDst) hn1AngDst->setVerboseLevel(verbose);
121  if (hn2AngDst) hn2AngDst->setVerboseLevel(verbose);
122  if (gnAngDst) gnAngDst->setVerboseLevel(verbose);
123  if (npAngDst) npAngDst->setVerboseLevel(verbose);
124  if (hn3BodyDst) hn3BodyDst->setVerboseLevel(verbose);
125  if (nn3BodyDst) nn3BodyDst->setVerboseLevel(verbose);
126 }
127 
128 
129 // Return appropriate distribution generator for specified interaction
130 
131 const G4VTwoBodyAngDst*
133  // TEMPORARY: Three-body distributions for hN/NN
134  if (fs==0 && kw==0) {
135  if (is == pro*pro || is == pro*neu || is == neu*neu) return nn3BodyDst;
136  else return hn3BodyDst;
137  }
138 
139  // gamma-nucleon -> nucleon pi0
140  if ((is == gam*pro && fs == pro*pi0) ||
141  (is == gam*neu && fs == neu*pi0)) {
142  return gp_ppi0;
143  }
144 
145  // gamma-nucleon charge exchange
146  if ((is == gam*pro && fs == neu*pip) ||
147  (is == gam*neu && fs == pro*pim)) {
148  return gp_npip;
149  }
150 
151  // pp and nn elastic
152  if (is == pro*pro || is == neu*neu) return ppAngDst;
153 
154  // np and pn elastic
155  if (is == pro*neu) return npAngDst;
156 
157  // pi+ p and pi- n elastic
158  if ((fs == is) && (is == pip*pro || is == pim*neu) ) return pippAngDst;
159 
160  // pi- p and pi+ n elastic
161  if ((fs == is) && (is == pim*pro || is == pip*neu) ) return pimpAngDst;
162 
163  // pi0 p and pi0 n elastic
164  if ((fs == is) && (is == pi0*pro || is == pi0*neu) ) return pi0pAngDst;
165 
166  // pi- p -> pi0 n, pi+ n -> pi0 p, pi0 p -> pi+ n, pi0 n -> pi- p
167  if ((is == pim*pro && fs == pi0*neu) || (is == pip*neu && fs == pi0*pip) ||
168  (is == pi0*pro && fs == pip*neu) || (is == pi0*neu && fs == pim*pro) )
169  return pipCXAngDst;
170 
171  // hyperon-nucleon
172  if (is == pro*lam || is == pro*sp || is == pro*s0 ||
173  is == pro*sm || is == pro*xi0 || is == pro*xim ||
174  is == pro*om ||
175  is == neu*lam || is == neu*sp || is == neu*s0 ||
176  is == neu*sm || is == neu*xi0 || is == neu*xim ||
177  is == neu*om) {
178  return nnAngDst;
179  }
180 
181  // gamma p -> K Y (and isospin variants)
182  if (kw == 2 && (is == pro*gam || is == neu*gam)) {
183  return gnAngDst;
184  }
185 
186  // pion-nucleon strangeness production
187  if (kw == 2) {
188  return qxAngDst;
189  }
190 
191  // gamma p, k+p, k0bp, gamma n, k-n, or k0n
192  if (is == pro*gam ||
193  is == pro*kpl || is == pro*k0b ||
194  is == neu*gam ||
195  is == neu*kmi || is == neu*k0) {
196  return hn1AngDst;
197  }
198 
199  // k-p, k0bn, k+n, or k0p
200  if (is == pro*kmi || is == pro*k0 ||
201  is == neu*kpl || is == neu*k0b) {
202  return hn2AngDst;
203  }
204 
205  // Invalid interaction
206  return 0;
207 }
G4NuclNucl3BodyAngDst * nn3BodyDst
G4PipP2PipPAngDst * pippAngDst
G4PimP2PimPAngDst * pimpAngDst
static G4ThreadLocal G4TwoBodyAngularDist * theInstance
G4GamP2PPi0AngDst * gp_ppi0
G4PiNInelasticAngDst * qxAngDst
static const G4TwoBodyAngularDist * GetInstance()
#define G4ThreadLocal
Definition: tls.hh:89
int G4int
Definition: G4Types.hh:78
G4HadNElastic2AngDst * hn2AngDst
G4NuclNuclAngDst * nnAngDst
G4GamP2NPipAngDst * gp_npip
void Register(T *inst)
Definition: G4AutoDelete.hh:65
const G4VTwoBodyAngDst * ChooseDist(G4int is, G4int fs, G4int kw) const
virtual void setVerboseLevel(G4int verbose=0)
G4Pi0P2Pi0PAngDst * pi0pAngDst
static void setVerboseLevel(G4int vb=0)
G4PimP2Pi0NAngDst * pipCXAngDst
G4HadNElastic1AngDst * hn1AngDst
G4GammaNuclAngDst * gnAngDst
G4HadNucl3BodyAngDst * hn3BodyDst
void passVerbose(G4int verbose)