Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4CascadeParamMessenger.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 // $Id: G4CascadeParamMessenger.cc 72016 2013-07-03 16:24:15Z mkelsey $
27 // Define simple UI commands as alternative to environment variables
28 //
29 // 20130304 M. Kelsey -- Add flag to collect and display cascade structure
30 // 20130621 M. Kelsey -- Add flag for CHECK_ECONS, replacing #ifdef's; add
31 // flag to use three-body momentum parametrizations
32 // 20130703 M. Kelsey -- Add flag for USE_PHASESPACE
33 // 20141030 M. Kelsey -- Add flag to enable direct pi-N absorption
34 // 20141211 M. Kelsey -- Change PIN_ABSORPTION flag to double, for energy cut
35 
37 #include "G4CascadeParameters.hh"
38 #include "G4UIcmdWithABool.hh"
39 #include "G4UIcmdWithADouble.hh"
40 #include "G4UIcmdWithAString.hh"
41 #include "G4UIcmdWithAnInteger.hh"
43 #include "G4UIcommand.hh"
44 #include "G4UIcommandTree.hh"
45 #include "G4UIdirectory.hh"
46 #include "G4UImanager.hh"
47 
48 
49 // Constructor and destructor
50 
52  : G4UImessenger(), theParams(params), cmdDir(0), localCmdDir(false) {
53  // NOTE: Put under same top-level tree as EM
54  CreateDirectory("/process/had/cascade/","Bertini-esque cascade parameters");
55 
56  verboseCmd = CreateCommand<G4UIcmdWithAnInteger>("verbose",
57  "Enable information messages");
58  balanceCmd = CreateCommand<G4UIcmdWithABool>("checkBalance",
59  "Enable internal conservation checking");
60  reportCmd = CreateCommand<G4UIcmdWithoutParameter>("report",
61  "Dump all non-default parameter settings");
62  usePreCoCmd = CreateCommand<G4UIcmdWithABool>("usePreCompound",
63  "Use PreCompoundModel for nuclear de-excitation");
64  doCoalCmd = CreateCommand<G4UIcmdWithABool>("doCoalescence",
65  "Apply final-state nucleon clustering");
66  piNAbsCmd = CreateCommand<G4UIcmdWithADouble>("piNAbsorption",
67  "Probability for pion absorption on single nucleon");
68  historyCmd = CreateCommand<G4UIcmdWithABool>("showHistory",
69  "Collect and report full structure of cascade");
70  use3BodyCmd = CreateCommand<G4UIcmdWithABool>("use3BodyMom",
71  "Use three-body momentum parametrizations");
72  usePSCmd = CreateCommand<G4UIcmdWithABool>("usePhaseSpace",
73  "Use Kopylov N-body momentum generator");
74  randomFileCmd = CreateCommand<G4UIcmdWithAString>("randomFile",
75  "Save random-engine to file at each interaction");
76  nucUseBestCmd = CreateCommand<G4UIcmdWithABool>("useBestNuclearModel",
77  "Use all physical-units for nuclear structure");
78  nucRad2parCmd = CreateCommand<G4UIcmdWithADouble>("useTwoParamNuclearRadius",
79  "Use R = C1*cbrt(A) + C2/cbrt(A)");
80  nucRadScaleCmd = CreateCommand<G4UIcmdWithADouble>("nuclearRadiusScale",
81  "Set length scale for nuclear model");
82  nucRadSmallCmd = CreateCommand<G4UIcmdWithADouble>("smallNucleusRadius",
83  "Set radius of A<4 nuclei");
84  nucRadAlphaCmd = CreateCommand<G4UIcmdWithADouble>("alphaRadiusScale",
85  "Fraction of small-radius for He-4");
86  nucRadTrailingCmd = CreateCommand<G4UIcmdWithADouble>("shadowningRadius",
87  "Effective nucleon radius for trailing effect");
88  nucFermiScaleCmd = CreateCommand<G4UIcmdWithADouble>("fermiScale",
89  "Scale factor for fermi momentum");
90  nucXsecScaleCmd = CreateCommand<G4UIcmdWithADouble>("crossSectionScale",
91  "Scale fator for total cross-sections");
92  nucGammaQDCmd = CreateCommand<G4UIcmdWithADouble>("gammaQuasiDeutScale",
93  "Scale factor for gamma-quasideutron cross-sections");
94  coalDPmax2Cmd = CreateCommand<G4UIcmdWithADouble>("cluster2DPmax",
95  "Maximum momentum for p-n clusters");
96  coalDPmax3Cmd = CreateCommand<G4UIcmdWithADouble>("cluster3DPmax",
97  "Maximum momentum for ppn/pnn clusters");
98  coalDPmax4Cmd = CreateCommand<G4UIcmdWithADouble>("cluster4DPmax",
99  "Maximum momentum for alpha clusters");
100 }
101 
103  delete verboseCmd;
104  delete balanceCmd;
105  delete reportCmd;
106  delete usePreCoCmd;
107  delete doCoalCmd;
108  delete piNAbsCmd;
109  delete historyCmd;
110  delete use3BodyCmd;
111  delete usePSCmd;
112  delete randomFileCmd;
113  delete nucUseBestCmd;
114  delete nucRad2parCmd;
115  delete nucRadScaleCmd;
116  delete nucRadSmallCmd;
117  delete nucRadAlphaCmd;
118  delete nucRadTrailingCmd;
119  delete nucFermiScaleCmd;
120  delete nucXsecScaleCmd;
121  delete nucGammaQDCmd;
122  delete coalDPmax2Cmd;
123  delete coalDPmax3Cmd;
124  delete coalDPmax4Cmd;
125  if (localCmdDir) delete cmdDir;
126 }
127 
128 
129 // Create or reuse existing UIdirectory path
130 
132  const char* desc) {
134  if (!UIman) return;
135 
136  // Directory path must be absolute, prepend "/" if ncessary
137  G4String fullPath = path;
138  if (fullPath(0) != '/') fullPath.prepend("/");
139  if (fullPath(fullPath.length()-1) != '/') fullPath.append("/");
140 
141  // See if input path has already been registered
142  G4UIcommand* foundPath = UIman->GetTree()->FindPath(fullPath);
143  if (foundPath) cmdDir = dynamic_cast<G4UIdirectory*>(foundPath);
144 
145  if (!cmdDir) { // Create local deletable directory
146  localCmdDir = true;
147  cmdDir = new G4UIdirectory(fullPath.c_str());
148  cmdDir->SetGuidance(desc);
149  }
150 }
151 
152 
153 // Use command argument (literal string) to set envvar maps in container
154 
156  if (cmd == reportCmd) theParams->DumpConfig(G4cout);
157 
158  if (cmd == verboseCmd)
159  theParams->G4CASCADE_VERBOSE = strdup(arg.c_str());
160 
161  if (cmd == balanceCmd)
162  theParams->G4CASCADE_CHECK_ECONS = StoB(arg) ? strdup(arg.c_str()) : 0;
163 
164  if (cmd == usePreCoCmd)
165  theParams->G4CASCADE_USE_PRECOMPOUND = StoB(arg) ? strdup(arg.c_str()) : 0;
166 
167  if (cmd == doCoalCmd)
168  theParams->G4CASCADE_DO_COALESCENCE = StoB(arg) ? strdup(arg.c_str()) : 0;
169 
170  if (cmd == piNAbsCmd)
171  theParams->G4CASCADE_PIN_ABSORPTION = strdup(arg.c_str());
172 
173  if (cmd == historyCmd)
174  theParams->G4CASCADE_SHOW_HISTORY = StoB(arg) ? strdup(arg.c_str()) : 0;
175 
176  if (cmd == use3BodyCmd)
177  theParams->G4CASCADE_USE_3BODYMOM = StoB(arg) ? strdup(arg.c_str()) : 0;
178 
179  if (cmd == usePSCmd)
180  theParams->G4CASCADE_USE_PHASESPACE = StoB(arg) ? strdup(arg.c_str()) : 0;
181 
182  if (cmd == randomFileCmd)
183  theParams->G4CASCADE_RANDOM_FILE = arg.empty() ? 0 : strdup(arg.c_str());
184 
185  if (cmd == nucUseBestCmd)
186  theParams->G4NUCMODEL_USE_BEST = StoB(arg) ? strdup(arg.c_str()) : 0;
187 
188 // if (cmd == nucRad2parCmd)
189 // theParams->G4NUCMODEL_RAD_2PAR = strdup(arg.c_str());
190  if (cmd == nucRad2parCmd)
191  theParams->G4NUCMODEL_RAD_2PAR = StoB(arg) ? strdup(arg.c_str()) : 0;
192 
193  if (cmd == nucRadScaleCmd)
194  theParams->G4NUCMODEL_RAD_SCALE = strdup(arg.c_str());
195 
196  if (cmd == nucRadSmallCmd)
197  theParams->G4NUCMODEL_RAD_SMALL = strdup(arg.c_str());
198 
199  if (cmd == nucRadAlphaCmd)
200  theParams->G4NUCMODEL_RAD_ALPHA = strdup(arg.c_str());
201 
202  if (cmd == nucRadTrailingCmd)
203  theParams->G4NUCMODEL_RAD_TRAILING = strdup(arg.c_str());
204 
205  if (cmd == nucFermiScaleCmd)
206  theParams->G4NUCMODEL_FERMI_SCALE = strdup(arg.c_str());
207 
208  if (cmd == nucXsecScaleCmd)
209  theParams->G4NUCMODEL_XSEC_SCALE = strdup(arg.c_str());
210 
211  if (cmd == nucGammaQDCmd)
212  theParams->G4NUCMODEL_GAMMAQD = strdup(arg.c_str());
213 
214  if (cmd == coalDPmax2Cmd)
215  theParams->DPMAX_2CLUSTER = strdup(arg.c_str());
216 
217  if (cmd == coalDPmax3Cmd)
218  theParams->DPMAX_3CLUSTER = strdup(arg.c_str());
219 
220  if (cmd == coalDPmax4Cmd)
221  theParams->DPMAX_4CLUSTER = strdup(arg.c_str());
222 
223  theParams->Initialize(); // Update numerical values from settings
224 }
const char * G4CASCADE_USE_PHASESPACE
const char * G4CASCADE_PIN_ABSORPTION
void SetGuidance(const char *aGuidance)
Definition: G4UIcommand.hh:161
G4UIcommandTree * GetTree() const
Definition: G4UImanager.hh:208
G4UIcmdWithADouble * nucRadSmallCmd
const char * G4NUCMODEL_RAD_SMALL
const char * G4CASCADE_USE_3BODYMOM
G4UIcmdWithADouble * nucFermiScaleCmd
virtual void SetNewValue(G4UIcommand *command, G4String newValue)
G4UIcmdWithADouble * coalDPmax4Cmd
G4UIcmdWithADouble * nucRad2parCmd
G4UIcmdWithADouble * nucGammaQDCmd
G4UIcmdWithADouble * coalDPmax3Cmd
const char * G4NUCMODEL_FERMI_SCALE
const char * G4NUCMODEL_GAMMAQD
static G4UImanager * GetUIpointer()
Definition: G4UImanager.cc:73
const char * G4NUCMODEL_XSEC_SCALE
const char * G4CASCADE_CHECK_ECONS
void DumpConfig(std::ostream &os) const
G4UIcmdWithADouble * coalDPmax2Cmd
G4UIcmdWithAnInteger * verboseCmd
const char * G4CASCADE_DO_COALESCENCE
const char * G4NUCMODEL_RAD_ALPHA
const char * G4NUCMODEL_RAD_TRAILING
G4CascadeParamMessenger(G4CascadeParameters *params)
G4UIcmdWithAString * randomFileCmd
G4UIcmdWithADouble * piNAbsCmd
void CreateDirectory(const char *path, const char *desc)
const char * G4NUCMODEL_RAD_SCALE
const char * G4CASCADE_RANDOM_FILE
const char * G4CASCADE_SHOW_HISTORY
G4UIcmdWithoutParameter * reportCmd
const char * G4NUCMODEL_USE_BEST
const char * G4CASCADE_USE_PRECOMPOUND
const char * G4CASCADE_VERBOSE
G4GLOB_DLL std::ostream G4cout
G4bool StoB(G4String s)
const char * G4NUCMODEL_RAD_2PAR
G4UIcmdWithADouble * nucRadAlphaCmd
G4String & prepend(const char *)
G4UIcmdWithADouble * nucXsecScaleCmd
G4String & append(const G4String &)
G4UIcommand * FindPath(const char *commandPath) const
G4CascadeParameters * theParams
G4UIcmdWithADouble * nucRadScaleCmd
G4UIcmdWithADouble * nucRadTrailingCmd