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