Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4CascadeParamMessenger Class Reference

#include <G4CascadeParamMessenger.hh>

Inheritance diagram for G4CascadeParamMessenger:
Collaboration diagram for G4CascadeParamMessenger:

Public Member Functions

 G4CascadeParamMessenger (G4CascadeParameters *params)
 
virtual ~G4CascadeParamMessenger ()
 
virtual void SetNewValue (G4UIcommand *command, G4String newValue)
 
- Public Member Functions inherited from G4UImessenger
 G4UImessenger ()
 
 G4UImessenger (const G4String &path, const G4String &dsc, G4bool commandsToBeBroadcasted=true)
 
virtual ~G4UImessenger ()
 
virtual G4String GetCurrentValue (G4UIcommand *command)
 
G4bool operator== (const G4UImessenger &messenger) const
 
G4bool CommandsShouldBeInMaster () const
 

Protected Member Functions

void CreateDirectory (const char *path, const char *desc)
 
template<class T >
T * CreateCommand (const G4String &cmd, const G4String &desc)
 
- Protected Member Functions inherited from G4UImessenger
G4String ItoS (G4int i)
 
G4String DtoS (G4double a)
 
G4String BtoS (G4bool b)
 
G4int StoI (G4String s)
 
G4double StoD (G4String s)
 
G4bool StoB (G4String s)
 
void AddUIcommand (G4UIcommand *newCommand)
 
void CreateDirectory (const G4String &path, const G4String &dsc, G4bool commandsToBeBroadcasted=true)
 
template<typename T >
T * CreateCommand (const G4String &cname, const G4String &dsc)
 

Additional Inherited Members

- Protected Attributes inherited from G4UImessenger
G4UIdirectorybaseDir
 
G4String baseDirName
 
G4bool commandsShouldBeInMaster
 

Detailed Description

Definition at line 52 of file G4CascadeParamMessenger.hh.

Constructor & Destructor Documentation

G4CascadeParamMessenger::G4CascadeParamMessenger ( G4CascadeParameters params)

Definition at line 51 of file G4CascadeParamMessenger.cc.

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 }
void CreateDirectory(const char *path, const char *desc)

Here is the call graph for this function:

G4CascadeParamMessenger::~G4CascadeParamMessenger ( )
virtual

Definition at line 102 of file G4CascadeParamMessenger.cc.

102  {
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 }

Member Function Documentation

template<class T >
T* G4CascadeParamMessenger::CreateCommand ( const G4String cmd,
const G4String desc 
)
protected
void G4CascadeParamMessenger::CreateDirectory ( const char *  path,
const char *  desc 
)
protected

Definition at line 131 of file G4CascadeParamMessenger.cc.

132  {
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 }
G4UIcommand * FindPath(const char *commandPath) const
static G4UImanager * GetUIpointer()
Definition: G4UImanager.cc:59
G4String & prepend(const char *)
void SetGuidance(const char *aGuidance)
Definition: G4UIcommand.hh:161
G4UIcommandTree * GetTree() const
Definition: G4UImanager.hh:206
G4String & append(const G4String &)

Here is the call graph for this function:

Here is the caller graph for this function:

void G4CascadeParamMessenger::SetNewValue ( G4UIcommand command,
G4String  newValue 
)
virtual

Reimplemented from G4UImessenger.

Definition at line 155 of file G4CascadeParamMessenger.cc.

155  {
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 
191  if (cmd == nucRadScaleCmd)
192  theParams->G4NUCMODEL_RAD_SCALE = strdup(arg.c_str());
193 
194  if (cmd == nucRadSmallCmd)
195  theParams->G4NUCMODEL_RAD_SMALL = strdup(arg.c_str());
196 
197  if (cmd == nucRadAlphaCmd)
198  theParams->G4NUCMODEL_RAD_ALPHA = strdup(arg.c_str());
199 
200  if (cmd == nucRadTrailingCmd)
201  theParams->G4NUCMODEL_RAD_TRAILING = strdup(arg.c_str());
202 
203  if (cmd == nucFermiScaleCmd)
204  theParams->G4NUCMODEL_FERMI_SCALE = strdup(arg.c_str());
205 
206  if (cmd == nucXsecScaleCmd)
207  theParams->G4NUCMODEL_XSEC_SCALE = strdup(arg.c_str());
208 
209  if (cmd == nucGammaQDCmd)
210  theParams->G4NUCMODEL_GAMMAQD = strdup(arg.c_str());
211 
212  if (cmd == coalDPmax2Cmd)
213  theParams->DPMAX_2CLUSTER = strdup(arg.c_str());
214 
215  if (cmd == coalDPmax3Cmd)
216  theParams->DPMAX_3CLUSTER = strdup(arg.c_str());
217 
218  if (cmd == coalDPmax4Cmd)
219  theParams->DPMAX_4CLUSTER = strdup(arg.c_str());
220 
221  theParams->Initialize(); // Update numerical values from settings
222 }
G4bool StoB(G4String s)
G4GLOB_DLL std::ostream G4cout

Here is the call graph for this function:


The documentation for this class was generated from the following files: