Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4CascadeParameters.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: G4CascadeParameters.cc 72016 2013-07-03 16:24:15Z mkelsey $
27 // Encapsulate all user-configurable parameters with associated envvars
28 //
29 // 20120912 M. Kelsey -- Add interface to support UI commands
30 // 20130304 M. Kelsey -- Add flag to collect and display cascade structure
31 // 20130308 M. Kelsey -- Add flag to use separate 3-body momentum generators
32 // 20130421 M. Kelsey -- Add flag for CHECK_ECONS, replacing #ifdef's
33 // 20130702 M. Kelsey -- Add flag to use N-body phase-space generator
34 // 20140311 G. Cosmo -- Implement standard (non-const) singleton pattern
35 // 20140929 M. Kelsey -- Enable some parameters as default true (must be set
36 // '0' for false): PreCompound, phase-space, clustering,
37 // trailing effect
38 // 20141111 M. Kelsey -- Revert defaults for PreCompound, phase-space,
39 // and trailing effect.
40 // 20141121 Use G4AutoDelete to avoid end-of-thread memory leaks
41 // 20141211 M. Kelsey -- Change PIN_ABSORPTION flag to double, for energy cut
42 
43 #include "G4CascadeParameters.hh"
45 #include "G4AutoDelete.hh"
46 #include <stdlib.h>
47 #include <iostream>
48 using std::endl;
49 
50 
51 // Singleton accessor
52 
53 G4CascadeParameters* G4CascadeParameters::fpInstance = 0;
54 
56  if (!fpInstance) {
57  fpInstance = new G4CascadeParameters;
58  G4AutoDelete::Register(fpInstance);
59  }
60 
61  return fpInstance;
62 }
63 
64 
65 // Constructor initializes everything once
66 
67 #define OLD_RADIUS_UNITS (3.3836/1.2) // Used with NucModel params
68 
69 G4CascadeParameters::G4CascadeParameters()
70  : G4CASCADE_VERBOSE(getenv("G4CASCADE_VERBOSE")),
71  G4CASCADE_CHECK_ECONS(getenv("G4CASCADE_CHECK_ECONS")),
72  G4CASCADE_USE_PRECOMPOUND(getenv("G4CASCADE_USE_PRECOMPOUND")),
73  G4CASCADE_DO_COALESCENCE(getenv("G4CASCADE_DO_COALESCENCE")),
74  G4CASCADE_SHOW_HISTORY(getenv("G4CASCADE_SHOW_HISTORY")),
75  G4CASCADE_USE_3BODYMOM(getenv("G4CASCADE_USE_3BODYMOM")),
76  G4CASCADE_USE_PHASESPACE(getenv("G4CASCADE_USE_PHASESPACE")),
77  G4CASCADE_PIN_ABSORPTION(getenv("G4CASCADE_PIN_ABSORPTION")),
78  G4CASCADE_RANDOM_FILE(getenv("G4CASCADE_RANDOM_FILE")),
79  G4NUCMODEL_USE_BEST(getenv("G4NUCMODEL_USE_BEST")),
80  G4NUCMODEL_RAD_2PAR(getenv("G4NUCMODEL_RAD_2PAR")),
81  G4NUCMODEL_RAD_SCALE(getenv("G4NUCMODEL_RAD_SCALE")),
82  G4NUCMODEL_RAD_SMALL(getenv("G4NUCMODEL_RAD_SMALL")),
83  G4NUCMODEL_RAD_ALPHA(getenv("G4NUCMODEL_RAD_ALPHA")),
84  G4NUCMODEL_RAD_TRAILING(getenv("G4NUCMODEL_RAD_TRAILING")),
85  G4NUCMODEL_FERMI_SCALE(getenv("G4NUCMODEL_FERMI_SCALE")),
86  G4NUCMODEL_XSEC_SCALE(getenv("G4NUCMODEL_XSEC_SCALE")),
87  G4NUCMODEL_GAMMAQD(getenv("G4NUCMODEL_GAMMAQD")),
88  DPMAX_2CLUSTER(getenv("DPMAX_2CLUSTER")),
89  DPMAX_3CLUSTER(getenv("DPMAX_3CLUSTER")),
90  DPMAX_4CLUSTER(getenv("DPMAX_4CLUSTER")),
91  messenger(0) {
92  messenger = new G4CascadeParamMessenger(this);
93  Initialize();
94 }
95 
96 void G4CascadeParameters::Initialize() {
97  VERBOSE_LEVEL = (G4CASCADE_VERBOSE ? atoi(G4CASCADE_VERBOSE) : 0);
98  CHECK_ECONS = (0!=G4CASCADE_CHECK_ECONS);
99  USE_PRECOMPOUND = (0!=G4CASCADE_USE_PRECOMPOUND &&
100  G4CASCADE_USE_PRECOMPOUND[0]!='0');
101  DO_COALESCENCE = (0==G4CASCADE_DO_COALESCENCE ||
102  G4CASCADE_DO_COALESCENCE[0]!='0');
103  SHOW_HISTORY = (0!=G4CASCADE_SHOW_HISTORY);
104  USE_3BODYMOM = (0!=G4CASCADE_USE_3BODYMOM);
105  USE_PHASESPACE = (0!=G4CASCADE_USE_PHASESPACE &&
106  G4CASCADE_USE_PHASESPACE[0]!='0');
107  PIN_ABSORPTION = (G4CASCADE_PIN_ABSORPTION ? strtod(G4CASCADE_PIN_ABSORPTION,0)
108  : 0.);
109  RANDOM_FILE = (G4CASCADE_RANDOM_FILE ? G4CASCADE_RANDOM_FILE : "");
110  BEST_PAR = (0!=G4NUCMODEL_USE_BEST);
111  TWOPARAM_RADIUS = (0!=G4NUCMODEL_RAD_2PAR);
112  RADIUS_SCALE = (G4NUCMODEL_RAD_SCALE ? strtod(G4NUCMODEL_RAD_SCALE,0)
113  : (BEST_PAR?1.0:OLD_RADIUS_UNITS));
114  RADIUS_SMALL = ((G4NUCMODEL_RAD_SMALL ? strtod(G4NUCMODEL_RAD_SMALL,0)
115  : (BEST_PAR?1.992:(8.0/OLD_RADIUS_UNITS))) * RADIUS_SCALE);
116  RADIUS_ALPHA = (G4NUCMODEL_RAD_ALPHA ? strtod(G4NUCMODEL_RAD_ALPHA,0)
117  : (BEST_PAR?0.84:0.70));
118  RADIUS_TRAILING = ((G4NUCMODEL_RAD_TRAILING ? strtod(G4NUCMODEL_RAD_TRAILING,0)
119  : 0.) * RADIUS_SCALE);
120  FERMI_SCALE = ((G4NUCMODEL_FERMI_SCALE ? strtod(G4NUCMODEL_FERMI_SCALE,0)
121  : (BEST_PAR?0.685:(1.932/OLD_RADIUS_UNITS))) * RADIUS_SCALE);
122  XSEC_SCALE = (G4NUCMODEL_XSEC_SCALE ? strtod(G4NUCMODEL_XSEC_SCALE,0)
123  : (BEST_PAR?0.1:1.0) );
124  GAMMAQD_SCALE = (G4NUCMODEL_GAMMAQD?strtod(G4NUCMODEL_GAMMAQD,0):1.);
125  DPMAX_DOUBLET = (DPMAX_2CLUSTER ? strtod(DPMAX_2CLUSTER,0) : 0.090);
126  DPMAX_TRIPLET = (DPMAX_3CLUSTER ? strtod(DPMAX_3CLUSTER,0) : 0.108);
127  DPMAX_ALPHA = (DPMAX_4CLUSTER ? strtod(DPMAX_4CLUSTER,0) : 0.115);
128 }
129 
131  delete messenger;
132 }
133 
134 
135 // Report any non-default parameters (used by G4CascadeInterface)
136 
137 void G4CascadeParameters::DumpConfig(std::ostream& os) const {
138  if (G4CASCADE_VERBOSE)
139  os << "G4CASCADE_VERBOSE = " << G4CASCADE_VERBOSE << endl;
140  if (G4CASCADE_CHECK_ECONS)
141  os << "G4CASCADE_CHECK_ECONS = " << G4CASCADE_CHECK_ECONS << endl;
142  if (G4CASCADE_USE_PRECOMPOUND)
143  os << "G4CASCADE_USE_PRECOMPOUND = " << G4CASCADE_USE_PRECOMPOUND << endl;
144  if (G4CASCADE_DO_COALESCENCE)
145  os << "G4CASCADE_DO_COALESCENCE = " << G4CASCADE_DO_COALESCENCE << endl;
146  if (G4CASCADE_PIN_ABSORPTION)
147  os << "G4CASCADE_PIN_ABSORPTION = " << G4CASCADE_PIN_ABSORPTION << endl;
148  if (G4CASCADE_SHOW_HISTORY)
149  os << "G4CASCADE_SHOW_HISTORY = " << G4CASCADE_SHOW_HISTORY << endl;
150  if (G4CASCADE_USE_3BODYMOM)
151  os << "G4CASCADE_USE_3BODYMOM = " << G4CASCADE_USE_3BODYMOM << endl;
152  if (G4CASCADE_USE_PHASESPACE)
153  os << "G4CASCADE_USE_PHASESPACE = " << G4CASCADE_USE_PHASESPACE << endl;
154  if (G4CASCADE_RANDOM_FILE)
155  os << "G4CASCADE_RANDOM_FILE = " << G4CASCADE_RANDOM_FILE << endl;
156  if (G4NUCMODEL_USE_BEST)
157  os << "G4NUCMODEL_USE_BEST = " << G4NUCMODEL_USE_BEST << endl;
158  if (G4NUCMODEL_RAD_2PAR)
159  os << "G4NUCMODEL_RAD_2PAR = " << G4NUCMODEL_RAD_2PAR << endl;
160  if (G4NUCMODEL_RAD_SCALE)
161  os << "G4NUCMODEL_RAD_SCALE = " << G4NUCMODEL_RAD_SCALE << endl;
162  if (G4NUCMODEL_RAD_SMALL)
163  os << "G4NUCMODEL_RAD_SMALL = " << G4NUCMODEL_RAD_SMALL << endl;
164  if (G4NUCMODEL_RAD_ALPHA)
165  os << "G4NUCMODEL_RAD_ALPHA = " << G4NUCMODEL_RAD_ALPHA << endl;
166  if (G4NUCMODEL_RAD_TRAILING)
167  os << "G4NUCMODEL_RAD_TRAILING = " << G4NUCMODEL_RAD_TRAILING << endl;
168  if (G4NUCMODEL_FERMI_SCALE)
169  os << "G4NUCMODEL_FERMI_SCALE = " << G4NUCMODEL_FERMI_SCALE << endl;
170  if (G4NUCMODEL_XSEC_SCALE)
171  os << "G4NUCMODEL_XSEC_SCALE = " << G4NUCMODEL_XSEC_SCALE << endl;
172  if (G4NUCMODEL_GAMMAQD)
173  os << "G4NUCMODEL_GAMMAQD = " << G4NUCMODEL_GAMMAQD << endl;
174  if (DPMAX_2CLUSTER)
175  os << "DPMAX_2CLUSTER = " << DPMAX_2CLUSTER << endl;
176  if (DPMAX_3CLUSTER)
177  os << "DPMAX_3CLUSTER = " << DPMAX_3CLUSTER << endl;
178  if (DPMAX_4CLUSTER)
179  os << "DPMAX_4CLUSTER = " << DPMAX_4CLUSTER << endl;
180 }
#define OLD_RADIUS_UNITS
void Register(T *inst)
Definition: G4AutoDelete.hh:65
static const G4CascadeParameters * Instance()