Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Protected Member Functions
G4CascadeParamMessenger Class Reference

#include <G4CascadeParamMessenger.hh>

Inheritance diagram for G4CascadeParamMessenger:
G4UImessenger

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
 

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
 

Detailed Description

Definition at line 50 of file G4CascadeParamMessenger.hh.

Constructor & Destructor Documentation

G4CascadeParamMessenger::G4CascadeParamMessenger ( G4CascadeParameters params)

Definition at line 49 of file G4CascadeParamMessenger.cc.

References CreateDirectory().

50  : G4UImessenger(), theParams(params), cmdDir(0), localCmdDir(false) {
51  // NOTE: Put under same top-level tree as EM
52  CreateDirectory("/process/had/cascade/","Bertini-esque cascade parameters");
53 
54  verboseCmd = CreateCommand<G4UIcmdWithAnInteger>("verbose",
55  "Enable information messages");
56  balanceCmd = CreateCommand<G4UIcmdWithABool>("checkBalance",
57  "Enable internal conservation checking");
58  reportCmd = CreateCommand<G4UIcmdWithoutParameter>("report",
59  "Dump all non-default parameter settings");
60  usePreCoCmd = CreateCommand<G4UIcmdWithABool>("usePreCompound",
61  "Use PreCompoundModel for nuclear de-excitation");
62  doCoalCmd = CreateCommand<G4UIcmdWithABool>("doCoalescence",
63  "Apply final-state nucleon clustering");
64  historyCmd = CreateCommand<G4UIcmdWithABool>("showHistory",
65  "Collect and report full structure of cascade");
66  use3BodyCmd = CreateCommand<G4UIcmdWithABool>("use3BodyMom",
67  "Use three-body momentum parametrizations");
68  usePSCmd = CreateCommand<G4UIcmdWithABool>("usePhaseSpace",
69  "Use Kopylov N-body momentum generator");
70  randomFileCmd = CreateCommand<G4UIcmdWithAString>("randomFile",
71  "Save random-engine to file at each interaction");
72  nucUseBestCmd = CreateCommand<G4UIcmdWithABool>("useBestNuclearModel",
73  "Use all physical-units for nuclear structure");
74  nucRad2parCmd = CreateCommand<G4UIcmdWithADouble>("useTwoParamNuclearRadius",
75  "Use R = C1*cbrt(A) + C2/cbrt(A)");
76  nucRadScaleCmd = CreateCommand<G4UIcmdWithADouble>("nuclearRadiusScale",
77  "Set length scale for nuclear model");
78  nucRadSmallCmd = CreateCommand<G4UIcmdWithADouble>("smallNucleusRadius",
79  "Set radius of A<4 nuclei");
80  nucRadAlphaCmd = CreateCommand<G4UIcmdWithADouble>("alphaRadiusScale",
81  "Fraction of small-radius for He-4");
82  nucRadTrailingCmd = CreateCommand<G4UIcmdWithADouble>("shadowningRadius",
83  "Effective nucleon radius for trailing effect");
84  nucFermiScaleCmd = CreateCommand<G4UIcmdWithADouble>("fermiScale",
85  "Scale factor for fermi momentum");
86  nucXsecScaleCmd = CreateCommand<G4UIcmdWithADouble>("crossSectionScale",
87  "Scale fator for total cross-sections");
88  nucGammaQDCmd = CreateCommand<G4UIcmdWithADouble>("gammaQuasiDeutScale",
89  "Scale factor for gamma-quasideutron cross-sections");
90  coalDPmax2Cmd = CreateCommand<G4UIcmdWithADouble>("cluster2DPmax",
91  "Maximum momentum for p-n clusters");
92  coalDPmax3Cmd = CreateCommand<G4UIcmdWithADouble>("cluster3DPmax",
93  "Maximum momentum for ppn/pnn clusters");
94  coalDPmax4Cmd = CreateCommand<G4UIcmdWithADouble>("cluster4DPmax",
95  "Maximum momentum for alpha clusters");
96 }
void CreateDirectory(const char *path, const char *desc)
G4CascadeParamMessenger::~G4CascadeParamMessenger ( )
virtual

Definition at line 98 of file G4CascadeParamMessenger.cc.

98  {
99  delete verboseCmd;
100  delete balanceCmd;
101  delete reportCmd;
102  delete usePreCoCmd;
103  delete doCoalCmd;
104  delete historyCmd;
105  delete use3BodyCmd;
106  delete usePSCmd;
107  delete randomFileCmd;
108  delete nucUseBestCmd;
109  delete nucRad2parCmd;
110  delete nucRadScaleCmd;
111  delete nucRadSmallCmd;
112  delete nucRadAlphaCmd;
113  delete nucRadTrailingCmd;
114  delete nucFermiScaleCmd;
115  delete nucXsecScaleCmd;
116  delete nucGammaQDCmd;
117  delete coalDPmax2Cmd;
118  delete coalDPmax3Cmd;
119  delete coalDPmax4Cmd;
120  if (localCmdDir) delete cmdDir;
121 }

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 126 of file G4CascadeParamMessenger.cc.

References G4String::append(), G4UIcommandTree::FindPath(), G4UImanager::GetTree(), G4UImanager::GetUIpointer(), G4String::prepend(), and G4UIcommand::SetGuidance().

Referenced by G4CascadeParamMessenger().

127  {
129  if (!UIman) return;
130 
131  // Directory path must be absolute, prepend "/" if ncessary
132  G4String fullPath = path;
133  if (fullPath(0) != '/') fullPath.prepend("/");
134  if (fullPath(fullPath.length()-1) != '/') fullPath.append("/");
135 
136  // See if input path has already been registered
137  G4UIcommand* foundPath = UIman->GetTree()->FindPath(fullPath);
138  if (foundPath) cmdDir = dynamic_cast<G4UIdirectory*>(foundPath);
139 
140  if (!cmdDir) { // Create local deletable directory
141  localCmdDir = true;
142  cmdDir = new G4UIdirectory(fullPath.c_str());
143  cmdDir->SetGuidance(desc);
144  }
145 }
G4UIcommand * FindPath(const char *commandPath) const
static G4UImanager * GetUIpointer()
Definition: G4UImanager.cc:58
G4String & prepend(const char *)
void SetGuidance(const char *aGuidance)
Definition: G4UIcommand.hh:161
G4UIcommandTree * GetTree() const
Definition: G4UImanager.hh:206
G4String & append(const G4String &)
void G4CascadeParamMessenger::SetNewValue ( G4UIcommand command,
G4String  newValue 
)
virtual

Reimplemented from G4UImessenger.

Definition at line 150 of file G4CascadeParamMessenger.cc.

References G4cout, and G4UImessenger::StoB().

150  {
151  if (cmd == reportCmd) theParams->DumpConfig(G4cout);
152 
153  if (cmd == verboseCmd)
154  theParams->G4CASCADE_VERBOSE = strdup(arg.c_str());
155 
156  if (cmd == balanceCmd)
157  theParams->G4CASCADE_CHECK_ECONS = StoB(arg) ? strdup(arg.c_str()) : 0;
158 
159  if (cmd == usePreCoCmd)
160  theParams->G4CASCADE_USE_PRECOMPOUND = StoB(arg) ? strdup(arg.c_str()) : 0;
161 
162  if (cmd == doCoalCmd)
163  theParams->G4CASCADE_DO_COALESCENCE = StoB(arg) ? strdup(arg.c_str()) : 0;
164 
165  if (cmd == historyCmd)
166  theParams->G4CASCADE_SHOW_HISTORY = StoB(arg) ? strdup(arg.c_str()) : 0;
167 
168  if (cmd == use3BodyCmd)
169  theParams->G4CASCADE_USE_3BODYMOM = StoB(arg) ? strdup(arg.c_str()) : 0;
170 
171  if (cmd == usePSCmd)
172  theParams->G4CASCADE_USE_PHASESPACE = StoB(arg) ? strdup(arg.c_str()) : 0;
173 
174  if (cmd == randomFileCmd)
175  theParams->G4CASCADE_RANDOM_FILE = arg.empty() ? 0 : strdup(arg.c_str());
176 
177  if (cmd == nucUseBestCmd)
178  theParams->G4NUCMODEL_USE_BEST = StoB(arg) ? strdup(arg.c_str()) : 0;
179 
180  if (cmd == nucRad2parCmd)
181  theParams->G4NUCMODEL_RAD_2PAR = strdup(arg.c_str());
182 
183  if (cmd == nucRadScaleCmd)
184  theParams->G4NUCMODEL_RAD_SCALE = strdup(arg.c_str());
185 
186  if (cmd == nucRadSmallCmd)
187  theParams->G4NUCMODEL_RAD_SMALL = strdup(arg.c_str());
188 
189  if (cmd == nucRadAlphaCmd)
190  theParams->G4NUCMODEL_RAD_ALPHA = strdup(arg.c_str());
191 
192  if (cmd == nucRadTrailingCmd)
193  theParams->G4NUCMODEL_RAD_TRAILING = strdup(arg.c_str());
194 
195  if (cmd == nucFermiScaleCmd)
196  theParams->G4NUCMODEL_FERMI_SCALE = strdup(arg.c_str());
197 
198  if (cmd == nucXsecScaleCmd)
199  theParams->G4NUCMODEL_XSEC_SCALE = strdup(arg.c_str());
200 
201  if (cmd == nucGammaQDCmd)
202  theParams->G4NUCMODEL_GAMMAQD = strdup(arg.c_str());
203 
204  if (cmd == coalDPmax2Cmd)
205  theParams->DPMAX_2CLUSTER = strdup(arg.c_str());
206 
207  if (cmd == coalDPmax3Cmd)
208  theParams->DPMAX_3CLUSTER = strdup(arg.c_str());
209 
210  if (cmd == coalDPmax4Cmd)
211  theParams->DPMAX_4CLUSTER = strdup(arg.c_str());
212 
213  theParams->Initialize(); // Update numerical values from settings
214 }
char cmd[1024]
Definition: tracer.cxx:25
G4bool StoB(G4String s)
G4GLOB_DLL std::ostream G4cout

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