Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4EmLowEParametersMessenger Class Reference

#include <G4EmLowEParametersMessenger.hh>

Inheritance diagram for G4EmLowEParametersMessenger:

Public Member Functions

 G4EmLowEParametersMessenger (G4EmLowEParameters *)
 ~G4EmLowEParametersMessenger () override
void SetNewValue (G4UIcommand *, G4String) override
G4EmLowEParametersMessengeroperator= (const G4EmLowEParametersMessenger &right)=delete
 G4EmLowEParametersMessenger (const G4EmLowEParametersMessenger &)=delete
Public Member Functions inherited from G4UImessenger
 G4UImessenger ()=default
 G4UImessenger (const G4String &path, const G4String &dsc, G4bool commandsToBeBroadcasted=true)
virtual ~G4UImessenger ()
virtual G4String GetCurrentValue (G4UIcommand *command)
G4bool CommandsShouldBeInMaster () const

Additional Inherited Members

Protected Member Functions inherited from G4UImessenger
G4String ItoS (G4int i)
G4String LtoS (G4long l)
G4String DtoS (G4double a)
G4String BtoS (G4bool b)
G4int StoI (const G4String &s)
G4long StoL (const G4String &s)
G4double StoD (const G4String &s)
G4bool StoB (const 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)
Protected Attributes inherited from G4UImessenger
G4UIdirectorybaseDir = nullptr
G4String baseDirName = ""
G4bool commandsShouldBeInMaster = false

Detailed Description

Definition at line 64 of file G4EmLowEParametersMessenger.hh.

Constructor & Destructor Documentation

◆ G4EmLowEParametersMessenger() [1/2]

G4EmLowEParametersMessenger::G4EmLowEParametersMessenger ( G4EmLowEParameters * ptr)
explicit

Definition at line 58 of file G4EmLowEParametersMessenger.cc.

59 : theParameters(ptr)
60{
61 deCmd = new G4UIcmdWithABool("/process/em/fluo",this);
62 deCmd->SetGuidance("Enable/disable atomic deexcitation");
63 deCmd->SetParameterName("fluoFlag",true);
64 deCmd->SetDefaultValue(false);
65 deCmd->AvailableForStates(G4State_PreInit,G4State_Init,G4State_Idle);
66 deCmd->SetToBeBroadcasted(false);
67
68 dirFluoCmd = new G4UIcmdWithABool("/process/em/fluoBearden",this);
69 dirFluoCmd->SetGuidance("Enable/disable usage of Bearden fluorescence files");
70 dirFluoCmd->SetParameterName("fluoBeardenFlag",true);
71 dirFluoCmd->SetDefaultValue(false);
72 dirFluoCmd->AvailableForStates(G4State_PreInit,G4State_Init);
73 dirFluoCmd->SetToBeBroadcasted(false);
74
75 dirFluoCmd1 = new G4UIcmdWithABool("/process/em/fluoANSTO",this);
76 dirFluoCmd1->SetGuidance("Enable/disable usage of ANSTO fluorescence files");
77 dirFluoCmd1->SetParameterName("fluoANSTOFlag",true);
78 dirFluoCmd1->SetDefaultValue(false);
79 dirFluoCmd1->AvailableForStates(G4State_PreInit,G4State_Init);
80 dirFluoCmd1->SetToBeBroadcasted(false);
81
82 auCmd = new G4UIcmdWithABool("/process/em/auger",this);
83 auCmd->SetGuidance("Enable/disable Auger electrons production");
84 auCmd->SetParameterName("augerFlag",true);
85 auCmd->SetDefaultValue(false);
86 auCmd->AvailableForStates(G4State_PreInit,G4State_Init,G4State_Idle);
87 auCmd->SetToBeBroadcasted(false);
88
89 auCascadeCmd = new G4UIcmdWithABool("/process/em/augerCascade",this);
90 auCascadeCmd->SetGuidance("Enable/disable simulation of cascade of Auger electrons");
91 auCascadeCmd->SetParameterName("augerCascadeFlag",true);
92 auCascadeCmd->SetDefaultValue(false);
93 auCascadeCmd->AvailableForStates(G4State_PreInit,G4State_Init,G4State_Idle);
94 auCascadeCmd->SetToBeBroadcasted(false);
95
96 pixeCmd = new G4UIcmdWithABool("/process/em/pixe",this);
97 pixeCmd->SetGuidance("Enable/disable PIXE simulation");
98 pixeCmd->SetParameterName("pixeFlag",true);
99 pixeCmd->SetDefaultValue(false);
100 pixeCmd->AvailableForStates(G4State_PreInit,G4State_Init,G4State_Idle);
101 pixeCmd->SetToBeBroadcasted(false);
102
103 dcutCmd = new G4UIcmdWithABool("/process/em/deexcitationIgnoreCut",this);
104 dcutCmd->SetGuidance("Enable/Disable usage of cuts in de-excitation module");
105 dcutCmd->SetParameterName("deexcut",true);
106 dcutCmd->SetDefaultValue(false);
107 dcutCmd->AvailableForStates(G4State_PreInit,G4State_Init,G4State_Idle);
108 dcutCmd->SetToBeBroadcasted(false);
109
110 dnafCmd = new G4UIcmdWithABool("/process/dna/UseDNAFast",this);
111 dnafCmd->SetGuidance("Enable usage of fast sampling for DNA models");
112 dnafCmd->SetParameterName("dnaf",true);
113 dnafCmd->SetDefaultValue(false);
114 dnafCmd->AvailableForStates(G4State_PreInit);
115 dnafCmd->SetToBeBroadcasted(false);
116
117 dnasCmd = new G4UIcmdWithABool("/process/dna/UseDNAStationary",this);
118 dnasCmd->SetGuidance("Enable usage of Stationary option for DNA models");
119 dnasCmd->SetParameterName("dnas",true);
120 dnasCmd->SetDefaultValue(false);
121 dnasCmd->AvailableForStates(G4State_PreInit);
122 dnasCmd->SetToBeBroadcasted(false);
123
124 dnamscCmd = new G4UIcmdWithABool("/process/dna/UseDNAElectronMsc",this);
125 dnamscCmd->SetGuidance("Enable usage of e- msc for DNA");
126 dnamscCmd->SetParameterName("dnamsc",true);
127 dnamscCmd->SetDefaultValue(false);
128 dnamscCmd->AvailableForStates(G4State_PreInit);
129 dnamscCmd->SetToBeBroadcasted(false);
130
131 direFluoCmd = new G4UIcmdWithAString("/process/em/fluoDirectory",this);
132 direFluoCmd->SetGuidance("The name of PIXE cross section");
133 direFluoCmd->SetParameterName("fluoDirectory",true);
134 direFluoCmd->SetCandidates("Default Bearden ANSTO XDB_EADL");
135 direFluoCmd->AvailableForStates(G4State_PreInit,G4State_Idle);
136 direFluoCmd->SetToBeBroadcasted(false);
137
138 pixeXsCmd = new G4UIcmdWithAString("/process/em/pixeXSmodel",this);
139 pixeXsCmd->SetGuidance("The name of PIXE cross section");
140 pixeXsCmd->SetParameterName("pixeXS",true);
141 pixeXsCmd->SetCandidates("ECPSSR_Analytical Empirical ECPSSR_FormFactor ECPSSR_ANSTO");
142 pixeXsCmd->AvailableForStates(G4State_PreInit,G4State_Idle);
143 pixeXsCmd->SetToBeBroadcasted(false);
144
145 pixeeXsCmd = new G4UIcmdWithAString("/process/em/pixeElecXSmodel",this);
146 pixeeXsCmd->SetGuidance("The name of PIXE cross section for electron");
147 pixeeXsCmd->SetParameterName("pixeEXS",true);
148 pixeeXsCmd->SetCandidates("ECPSSR_Analytical Empirical Livermore Penelope");
149 pixeeXsCmd->AvailableForStates(G4State_PreInit,G4State_Idle);
150 pixeeXsCmd->SetToBeBroadcasted(false);
151
152 livCmd = new G4UIcmdWithAString("/process/em/LivermoreData",this);
153 livCmd->SetGuidance("The name of Livermore data directory");
154 livCmd->SetParameterName("livDir",true);
155 livCmd->SetCandidates("livermore epics_2017");
156 livCmd->AvailableForStates(G4State_PreInit);
157 livCmd->SetToBeBroadcasted(false);
158
159 dnaSolCmd = new G4UIcmdWithAString("/process/dna/e-SolvationSubType",this);
160 dnaSolCmd->SetGuidance("The name of e- solvation DNA model");
161 dnaSolCmd->SetParameterName("dnaSol",true);
162 dnaSolCmd->SetCandidates("Ritchie1994 Terrisol1990 Meesungnoen2002 Kreipl2009 Meesungnoen2002_amorphous");
163 dnaSolCmd->AvailableForStates(G4State_PreInit);
164 dnaSolCmd->SetToBeBroadcasted(false);
165
166 dnaChemModel = new G4UIcmdWithAString("/process/chem/TimeStepModel",this);
167 dnaChemModel->SetGuidance("The name of DNA chemistry time step model");
168 dnaChemModel->SetParameterName("TimeStepModel",true);
169 dnaChemModel->SetCandidates("SBS IRT IRT_syn");
170 dnaChemModel->AvailableForStates(G4State_PreInit);
171 dnaChemModel->SetToBeBroadcasted(false);
172
173 meCmd = new G4UIcmdWithAString("/process/em/AddMicroElecRegion",this);
174 meCmd->SetGuidance("Activate MicroElec model in the G4Region");
175 meCmd->SetParameterName("MicroElec",true);
176 meCmd->AvailableForStates(G4State_PreInit);
177 meCmd->SetToBeBroadcasted(false);
178
179 dnaCmd = new G4UIcommand("/process/em/AddDNARegion",this);
180 dnaCmd->SetGuidance("Activate DNA in a G4Region.");
181 dnaCmd->SetGuidance(" regName : G4Region name");
182 dnaCmd->SetGuidance(" dnaType : DNA_opt0, DNA_Opt2, DNA_Opt4, DNA_Opt4a, DNA_Opt6, DNA_Opt6a, DNA_Opt7");
183 dnaCmd->AvailableForStates(G4State_PreInit);
184 dnaCmd->SetToBeBroadcasted(false);
185
186 auto regName = new G4UIparameter("regName",'s',false);
187 dnaCmd->SetParameter(regName);
188
189 auto type = new G4UIparameter("dnaType",'s',false);
190 dnaCmd->SetParameter(type);
191 type->SetParameterCandidates("DNA_Opt0 DNA_Opt2 DNA_Opt4 DNA_Opt4a DNA_Opt6 DNA_Opt6a DNA_Opt7");
192
193 deexCmd = new G4UIcommand("/process/em/deexcitation",this);
194 deexCmd->SetGuidance("Set deexcitation flags per G4Region.");
195 deexCmd->SetGuidance(" regName : G4Region name");
196 deexCmd->SetGuidance(" flagFluo : Fluorescence");
197 deexCmd->SetGuidance(" flagAuger : Auger");
198 deexCmd->SetGuidance(" flagPIXE : PIXE");
199 deexCmd->AvailableForStates(G4State_PreInit,G4State_Init,G4State_Idle);
200 deexCmd->SetToBeBroadcasted(false);
201
202 auto regNameD = new G4UIparameter("regName",'s',false);
203 deexCmd->SetParameter(regNameD);
204
205 auto flagFluo = new G4UIparameter("flagFluo",'s',false);
206 deexCmd->SetParameter(flagFluo);
207
208 auto flagAuger = new G4UIparameter("flagAuger",'s',false);
209 deexCmd->SetParameter(flagAuger);
210
211 auto flagPIXE = new G4UIparameter("flagPIXE",'s',false);
212 deexCmd->SetParameter(flagPIXE);
213
214}
@ G4State_Init
@ G4State_Idle
@ G4State_PreInit

Referenced by G4EmLowEParametersMessenger(), and operator=().

◆ ~G4EmLowEParametersMessenger()

G4EmLowEParametersMessenger::~G4EmLowEParametersMessenger ( )
override

Definition at line 218 of file G4EmLowEParametersMessenger.cc.

219{
220 delete deCmd;
221 delete dirFluoCmd;
222 delete dirFluoCmd1;
223 delete auCmd;
224 delete auCascadeCmd;
225 delete pixeCmd;
226 delete dcutCmd;
227 delete dnafCmd;
228 delete dnasCmd;
229 delete dnamscCmd;
230 delete pixeXsCmd;
231 delete pixeeXsCmd;
232 delete livCmd;
233 delete dnaSolCmd;
234 delete dnaChemModel;
235 delete direFluoCmd;
236 delete meCmd;
237 delete dnaCmd;
238 delete deexCmd;
239}

◆ G4EmLowEParametersMessenger() [2/2]

G4EmLowEParametersMessenger::G4EmLowEParametersMessenger ( const G4EmLowEParametersMessenger & )
delete

Member Function Documentation

◆ operator=()

G4EmLowEParametersMessenger & G4EmLowEParametersMessenger::operator= ( const G4EmLowEParametersMessenger & right)
delete

◆ SetNewValue()

void G4EmLowEParametersMessenger::SetNewValue ( G4UIcommand * command,
G4String newValue )
overridevirtual

Reimplemented from G4UImessenger.

Definition at line 243 of file G4EmLowEParametersMessenger.cc.

245{
246 G4bool physicsModified = false;
247 if (command == deCmd) {
248 theParameters->SetFluo(deCmd->GetNewBoolValue(newValue));
249 physicsModified = true;
250 } else if (command == dirFluoCmd) {
251 theParameters->SetBeardenFluoDir(dirFluoCmd->GetNewBoolValue(newValue));
252 physicsModified = true;
253 } else if (command == dirFluoCmd1) {
254 theParameters->SetANSTOFluoDir(dirFluoCmd1->GetNewBoolValue(newValue));
255 physicsModified = true;
256 } else if (command == auCmd) {
257 theParameters->SetAuger(auCmd->GetNewBoolValue(newValue));
258 physicsModified = true;
259 } else if (command == auCascadeCmd) {
260 theParameters->SetAuger(auCascadeCmd->GetNewBoolValue(newValue));
261 physicsModified = true;
262 } else if (command == pixeCmd) {
263 theParameters->SetPixe(pixeCmd->GetNewBoolValue(newValue));
264 physicsModified = true;
265 } else if (command == dcutCmd) {
266 theParameters->SetDeexcitationIgnoreCut(dcutCmd->GetNewBoolValue(newValue));
267 physicsModified = true;
268 } else if (command == dnafCmd) {
269 theParameters->SetDNAFast(dnafCmd->GetNewBoolValue(newValue));
270 } else if (command == dnasCmd) {
271 theParameters->SetDNAStationary(dnasCmd->GetNewBoolValue(newValue));
272 } else if (command == dnamscCmd) {
273 theParameters->SetDNAElectronMsc(dnamscCmd->GetNewBoolValue(newValue));
274 } else if (command == dnaSolCmd) {
276 if(newValue == "Ritchie1994") {
278 } else if(newValue == "Terrisol1990") {
280 } else if (newValue == "Meesungnoen2002") {
282 } else if (newValue == "Meesungnoen2002_amorphous") {
284 } else if (newValue == "Kreipl2009") {
286 }
287 theParameters->SetDNAeSolvationSubType(ttt);
288 } else if (command == dnaChemModel) {
290 if(newValue == "IRT") {
292 } else if(newValue == "SBS") {
294 } else if (newValue == "IRT_syn") {
296 }
297 theParameters->SetChemTimeStepModel(stepM);
298 } else if (command == direFluoCmd) {
300 if(newValue == "Bearden") { ttt = fluoBearden; }
301 else if(newValue == "ANSTO") { ttt = fluoANSTO; }
302 else if(newValue == "XDB_EADL") { ttt = fluoXDB_EADL; }
303 theParameters->SetFluoDirectory(ttt);
304 } else if (command == pixeXsCmd) {
305 theParameters->SetPIXECrossSectionModel(newValue);
306 physicsModified = true;
307 } else if (command == pixeeXsCmd) {
308 theParameters->SetPIXEElectronCrossSectionModel(newValue);
309 physicsModified = true;
310 } else if (command == livCmd) {
311 theParameters->SetLivermoreDataDir(newValue);
312 } else if (command == meCmd) {
313 theParameters->AddMicroElec(newValue);
314 } else if (command == dnaCmd) {
315 G4String s1(""),s2("");
316 std::istringstream is(newValue);
317 is >> s1 >> s2;
318 theParameters->AddDNA(s1, s2);
319 } else if (command == deexCmd) {
320 G4String s1 (""), s2(""), s3(""), s4("");
321 G4bool b2(false), b3(false), b4(false);
322 std::istringstream is(newValue);
323 is >> s1 >> s2 >> s3 >> s4;
324 if(s2 == "true") { b2 = true; }
325 if(s3 == "true") { b3 = true; }
326 if(s4 == "true") { b4 = true; }
327 theParameters->SetDeexActiveRegion(s1,b2,b3,b4);
328 physicsModified = true;
329 }
330
331 if(physicsModified) {
332 G4UImanager::GetUIpointer()->ApplyCommand("/run/physicsModified");
333 }
334}
G4ChemTimeStepModel
G4DNAModelSubType
@ fMeesungnoen2002eSolvation
@ fKreipl2009eSolvation
@ fDNAUnknownModel
@ fMeesungnoensolid2002eSolvation
@ fRitchie1994eSolvation
@ fTerrisol1990eSolvation
G4EmFluoDirectory
@ fluoBearden
@ fluoXDB_EADL
@ fluoDefault
@ fluoANSTO
bool G4bool
Definition G4Types.hh:86
G4int ApplyCommand(const char *aCommand)
static G4UImanager * GetUIpointer()

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