Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups 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 
35 #include "G4CascadeParameters.hh"
37 #include <stdlib.h>
38 #include <iostream>
39 using std::endl;
40 
41 
42 // Singleton accessor
43 
45  static const G4CascadeParameters theInstance;
46  return &theInstance;
47 }
48 
49 
50 // Constructor initializes everything once
51 
52 #define OLD_RADIUS_UNITS (3.3836/1.2) // Used with NucModel params
53 
54 G4CascadeParameters::G4CascadeParameters()
55  : G4CASCADE_VERBOSE(getenv("G4CASCADE_VERBOSE")),
56  G4CASCADE_CHECK_ECONS(getenv("G4CASCADE_CHECK_ECONS")),
57  G4CASCADE_USE_PRECOMPOUND(getenv("G4CASCADE_USE_PRECOMPOUND")),
58  G4CASCADE_DO_COALESCENCE(getenv("G4CASCADE_DO_COALESCENCE")),
59  G4CASCADE_SHOW_HISTORY(getenv("G4CASCADE_SHOW_HISTORY")),
60  G4CASCADE_USE_3BODYMOM(getenv("G4CASCADE_USE_3BODYMOM")),
61  G4CASCADE_USE_PHASESPACE(getenv("G4CASCADE_USE_PHASESPACE")),
62  G4CASCADE_RANDOM_FILE(getenv("G4CASCADE_RANDOM_FILE")),
63  G4NUCMODEL_USE_BEST(getenv("G4NUCMODEL_USE_BEST")),
64  G4NUCMODEL_RAD_2PAR(getenv("G4NUCMODEL_RAD_2PAR")),
65  G4NUCMODEL_RAD_SCALE(getenv("G4NUCMODEL_RAD_SCALE")),
66  G4NUCMODEL_RAD_SMALL(getenv("G4NUCMODEL_RAD_SMALL")),
67  G4NUCMODEL_RAD_ALPHA(getenv("G4NUCMODEL_RAD_ALPHA")),
68  G4NUCMODEL_RAD_TRAILING(getenv("G4NUCMODEL_RAD_TRAILING")),
69  G4NUCMODEL_FERMI_SCALE(getenv("G4NUCMODEL_FERMI_SCALE")),
70  G4NUCMODEL_XSEC_SCALE(getenv("G4NUCMODEL_XSEC_SCALE")),
71  G4NUCMODEL_GAMMAQD(getenv("G4NUCMODEL_GAMMAQD")),
72  DPMAX_2CLUSTER(getenv("DPMAX_2CLUSTER")),
73  DPMAX_3CLUSTER(getenv("DPMAX_3CLUSTER")),
74  DPMAX_4CLUSTER(getenv("DPMAX_4CLUSTER")),
75  messenger(0) {
76  messenger = new G4CascadeParamMessenger(this);
77  Initialize();
78 }
79 
80 void G4CascadeParameters::Initialize() {
81  VERBOSE_LEVEL = (G4CASCADE_VERBOSE ? atoi(G4CASCADE_VERBOSE) : 0);
82  CHECK_ECONS = (0!=G4CASCADE_CHECK_ECONS);
83  USE_PRECOMPOUND = (0!=G4CASCADE_USE_PRECOMPOUND);
84  DO_COALESCENCE = (0!=G4CASCADE_DO_COALESCENCE);
85  SHOW_HISTORY = (0!=G4CASCADE_SHOW_HISTORY);
86  USE_3BODYMOM = (0!=G4CASCADE_USE_3BODYMOM);
87  USE_PHASESPACE = (0!=G4CASCADE_USE_PHASESPACE);
88  RANDOM_FILE = (G4CASCADE_RANDOM_FILE ? G4CASCADE_RANDOM_FILE : "");
89  BEST_PAR = (0!=G4NUCMODEL_USE_BEST);
90  TWOPARAM_RADIUS = (0!=G4NUCMODEL_RAD_2PAR);
91  RADIUS_SCALE = (G4NUCMODEL_RAD_SCALE ? strtod(G4NUCMODEL_RAD_SCALE,0)
92  : (BEST_PAR?1.0:OLD_RADIUS_UNITS));
93  RADIUS_SMALL = ((G4NUCMODEL_RAD_SMALL ? strtod(G4NUCMODEL_RAD_SMALL,0)
94  : (BEST_PAR?1.992:(8.0/OLD_RADIUS_UNITS))) * RADIUS_SCALE);
95  RADIUS_ALPHA = (G4NUCMODEL_RAD_ALPHA ? strtod(G4NUCMODEL_RAD_ALPHA,0)
96  : (BEST_PAR?0.84:0.70));
97  RADIUS_TRAILING = ((G4NUCMODEL_RAD_TRAILING ? strtod(G4NUCMODEL_RAD_TRAILING,0)
98  : (BEST_PAR?0.70:0.0)) * RADIUS_SCALE);
99  FERMI_SCALE = ((G4NUCMODEL_FERMI_SCALE ? strtod(G4NUCMODEL_FERMI_SCALE,0)
100  : (BEST_PAR?0.685:(1.932/OLD_RADIUS_UNITS))) * RADIUS_SCALE);
101  XSEC_SCALE = (G4NUCMODEL_XSEC_SCALE ? strtod(G4NUCMODEL_XSEC_SCALE,0)
102  : (BEST_PAR?0.1:1.0) );
103  GAMMAQD_SCALE = (G4NUCMODEL_GAMMAQD?strtod(G4NUCMODEL_GAMMAQD,0):1.);
104  DPMAX_DOUBLET = (DPMAX_2CLUSTER ? strtod(DPMAX_2CLUSTER,0) : 0.090);
105  DPMAX_TRIPLET = (DPMAX_3CLUSTER ? strtod(DPMAX_3CLUSTER,0) : 0.108);
106  DPMAX_ALPHA = (DPMAX_4CLUSTER ? strtod(DPMAX_4CLUSTER,0) : 0.115);
107 }
108 
110  delete messenger;
111 }
112 
113 
114 // Report any non-default parameters (used by G4CascadeInterface)
115 
116 void G4CascadeParameters::DumpConfig(std::ostream& os) const {
117  if (G4CASCADE_VERBOSE)
118  os << "G4CASCADE_VERBOSE = " << G4CASCADE_VERBOSE << endl;
119  if (G4CASCADE_CHECK_ECONS)
120  os << "G4CASCADE_CHECK_ECONS = " << G4CASCADE_CHECK_ECONS << endl;
121  if (G4CASCADE_USE_PRECOMPOUND)
122  os << "G4CASCADE_USE_PRECOMPOUND = " << G4CASCADE_USE_PRECOMPOUND << endl;
123  if (G4CASCADE_DO_COALESCENCE)
124  os << "G4CASCADE_DO_COALESCENCE = " << G4CASCADE_DO_COALESCENCE << endl;
125  if (G4CASCADE_SHOW_HISTORY)
126  os << "G4CASCADE_SHOW_HISTORY = " << G4CASCADE_SHOW_HISTORY << endl;
127  if (G4CASCADE_USE_3BODYMOM)
128  os << "G4CASCADE_USE_3BODYMOM = " << G4CASCADE_USE_3BODYMOM << endl;
129  if (G4CASCADE_USE_PHASESPACE)
130  os << "G4CASCADE_USE_PHASESPACE = " << G4CASCADE_USE_PHASESPACE << endl;
131  if (G4CASCADE_RANDOM_FILE)
132  os << "G4CASCADE_RANDOM_FILE = " << G4CASCADE_RANDOM_FILE << endl;
133  if (G4NUCMODEL_USE_BEST)
134  os << "G4NUCMODEL_USE_BEST = " << G4NUCMODEL_USE_BEST << endl;
135  if (G4NUCMODEL_RAD_2PAR)
136  os << "G4NUCMODEL_RAD_2PAR = " << G4NUCMODEL_RAD_2PAR << endl;
137  if (G4NUCMODEL_RAD_SCALE)
138  os << "G4NUCMODEL_RAD_SCALE = " << G4NUCMODEL_RAD_SCALE << endl;
139  if (G4NUCMODEL_RAD_SMALL)
140  os << "G4NUCMODEL_RAD_SMALL = " << G4NUCMODEL_RAD_SMALL << endl;
141  if (G4NUCMODEL_RAD_ALPHA)
142  os << "G4NUCMODEL_RAD_ALPHA = " << G4NUCMODEL_RAD_ALPHA << endl;
143  if (G4NUCMODEL_RAD_TRAILING)
144  os << "G4NUCMODEL_RAD_TRAILING = " << G4NUCMODEL_RAD_TRAILING << endl;
145  if (G4NUCMODEL_FERMI_SCALE)
146  os << "G4NUCMODEL_FERMI_SCALE = " << G4NUCMODEL_FERMI_SCALE << endl;
147  if (G4NUCMODEL_XSEC_SCALE)
148  os << "G4NUCMODEL_XSEC_SCALE = " << G4NUCMODEL_XSEC_SCALE << endl;
149  if (G4NUCMODEL_GAMMAQD)
150  os << "G4NUCMODEL_GAMMAQD = " << G4NUCMODEL_GAMMAQD << endl;
151  if (DPMAX_2CLUSTER)
152  os << "DPMAX_2CLUSTER = " << DPMAX_2CLUSTER << endl;
153  if (DPMAX_3CLUSTER)
154  os << "DPMAX_3CLUSTER = " << DPMAX_3CLUSTER << endl;
155  if (DPMAX_4CLUSTER)
156  os << "DPMAX_4CLUSTER = " << DPMAX_4CLUSTER << endl;
157 }
#define OLD_RADIUS_UNITS
void Initialize()
Definition: errprop.cc:96
static const G4CascadeParameters * Instance()