Geant4-11
Public Member Functions | Protected Member Functions | Protected Attributes | Private Attributes
G4EmExtraParametersMessenger Class Reference

#include <G4EmExtraParametersMessenger.hh>

Inheritance diagram for G4EmExtraParametersMessenger:
G4UImessenger

Public Member Functions

G4bool CommandsShouldBeInMaster () const
 
 G4EmExtraParametersMessenger (const G4EmExtraParametersMessenger &)=delete
 
 G4EmExtraParametersMessenger (G4EmExtraParameters *)
 
virtual G4String GetCurrentValue (G4UIcommand *command)
 
G4bool operator!= (const G4UImessenger &messenger) const
 
G4EmExtraParametersMessengeroperator= (const G4EmExtraParametersMessenger &right)=delete
 
G4bool operator== (const G4UImessenger &messenger) const
 
void SetNewValue (G4UIcommand *, G4String) override
 
 ~G4EmExtraParametersMessenger () override
 

Protected Member Functions

void AddUIcommand (G4UIcommand *newCommand)
 
G4String BtoS (G4bool b)
 
template<typename T >
T * CreateCommand (const G4String &cname, const G4String &dsc)
 
void CreateDirectory (const G4String &path, const G4String &dsc, G4bool commandsToBeBroadcasted=true)
 
G4String DtoS (G4double a)
 
G4String ItoS (G4int i)
 
G4bool StoB (G4String s)
 
G4double StoD (G4String s)
 
G4int StoI (G4String s)
 
G4long StoL (G4String s)
 

Protected Attributes

G4UIdirectorybaseDir = nullptr
 
G4String baseDirName = ""
 
G4bool commandsShouldBeInMaster = false
 

Private Attributes

G4UIcommandbfCmd
 
G4UIcommandbsCmd
 
G4UIcmdWithABooldirSplitCmd
 
G4UIcmdWithADoubleAndUnitdirSplitRadiusCmd
 
G4UIcmdWith3VectorAndUnitdirSplitTargetCmd
 
G4UIcommandfiCmd
 
G4UIcommandmscoCmd
 
G4UIcommandpaiCmd
 
G4UIcmdWithABoolqeCmd
 
G4UIcommandStepFuncCmd
 
G4UIcommandStepFuncCmd1
 
G4UIcommandStepFuncCmd2
 
G4UIcommandStepFuncCmd3
 
G4UIcmdWithAStringSubSecCmd
 
G4EmExtraParameterstheParameters
 

Detailed Description

Definition at line 64 of file G4EmExtraParametersMessenger.hh.

Constructor & Destructor Documentation

◆ G4EmExtraParametersMessenger() [1/2]

G4EmExtraParametersMessenger::G4EmExtraParametersMessenger ( G4EmExtraParameters ptr)
explicit

Definition at line 58 of file G4EmExtraParametersMessenger.cc.

59 : theParameters(ptr)
60{
61 paiCmd = new G4UIcommand("/process/em/AddPAIRegion",this);
62 paiCmd->SetGuidance("Activate PAI in the G4Region.");
63 paiCmd->SetGuidance(" partName : particle name (default - all)");
64 paiCmd->SetGuidance(" regName : G4Region name");
65 paiCmd->SetGuidance(" paiType : PAI, PAIphoton");
68
69 G4UIparameter* part = new G4UIparameter("partName",'s',false);
70 paiCmd->SetParameter(part);
71
72 G4UIparameter* pregName = new G4UIparameter("regName",'s',false);
73 paiCmd->SetParameter(pregName);
74
75 G4UIparameter* ptype = new G4UIparameter("type",'s',false);
76 paiCmd->SetParameter(ptype);
77 ptype->SetParameterCandidates("pai PAI PAIphoton");
78
79 mscoCmd = new G4UIcommand("/process/em/AddEmRegion",this);
80 mscoCmd->SetGuidance("Add optional EM configuration for a G4Region.");
81 mscoCmd->SetGuidance(" regName : G4Region name");
82 mscoCmd->SetGuidance(" emType : G4EmStandard, G4EmStandard_opt1, ...");
84
85 G4UIparameter* mregName = new G4UIparameter("regName",'s',false);
86 mscoCmd->SetParameter(mregName);
87
88 G4UIparameter* mtype = new G4UIparameter("mscType",'s',false);
89 mscoCmd->SetParameter(mtype);
90 mtype->SetParameterCandidates("G4EmStandard G4EmStandard_opt1 G4EmStandard_opt2 G4EmStandard_opt3 G4EmStandard_opt4 G4EmStandardGS G4EmStandardSS G4EmLivermore G4EmPenelope G4RadioactiveDecay");
91
92 SubSecCmd = new G4UIcmdWithAString("/process/eLoss/subsecRegion",this);
93 SubSecCmd->SetGuidance("Enable subcut generation per region.");
94 SubSecCmd->SetGuidance(" Region : region name");
97
98 StepFuncCmd = new G4UIcommand("/process/eLoss/StepFunction",this);
99 StepFuncCmd->SetGuidance("Set the energy loss step limitation parameters for e+-.");
100 StepFuncCmd->SetGuidance(" dRoverR : max Range variation per step");
101 StepFuncCmd->SetGuidance(" finalRange: range for final step");
102 StepFuncCmd->SetGuidance(" unit : unit of finalRange");
105
106 G4UIparameter* dRoverRPrm = new G4UIparameter("dRoverR",'d',false);
107 dRoverRPrm->SetParameterRange("dRoverR>0. && dRoverR<=1.");
108 StepFuncCmd->SetParameter(dRoverRPrm);
109
110 G4UIparameter* finalRangePrm = new G4UIparameter("finalRange",'d',false);
111 finalRangePrm->SetParameterRange("finalRange>0.");
112 StepFuncCmd->SetParameter(finalRangePrm);
113
114 G4UIparameter* unitPrm = new G4UIparameter("unit",'s',true);
115 unitPrm->SetDefaultUnit("mm");
116 StepFuncCmd->SetParameter(unitPrm);
117
118 StepFuncCmd1 = new G4UIcommand("/process/eLoss/StepFunctionMuHad",this);
119 StepFuncCmd1->SetGuidance("Set the energy loss step limitation parameters for muon/hadron.");
120 StepFuncCmd1->SetGuidance(" dRoverR : max Range variation per step");
121 StepFuncCmd1->SetGuidance(" finalRange: range for final step");
124
125 G4UIparameter* dRoverRPrm1 = new G4UIparameter("dRoverRMuHad",'d',false);
126 dRoverRPrm1->SetParameterRange("dRoverRMuHad>0. && dRoverRMuHad<=1.");
127 StepFuncCmd1->SetParameter(dRoverRPrm1);
128
129 G4UIparameter* finalRangePrm1 = new G4UIparameter("finalRangeMuHad",'d',false);
130 finalRangePrm1->SetParameterRange("finalRangeMuHad>0.");
131 StepFuncCmd1->SetParameter(finalRangePrm1);
132
133 G4UIparameter* unitPrm1 = new G4UIparameter("unit",'s',true);
134 unitPrm1->SetDefaultValue("mm");
135 StepFuncCmd1->SetParameter(unitPrm1);
136
137 StepFuncCmd2 = new G4UIcommand("/process/eLoss/StepFunctionLightIons",this);
138 StepFuncCmd2->SetGuidance("Set the energy loss step limitation parameters for light ions.");
139 StepFuncCmd2->SetGuidance(" dRoverR : max Range variation per step");
140 StepFuncCmd2->SetGuidance(" finalRange: range for final step");
143
144 G4UIparameter* dRoverRPrm2 = new G4UIparameter("dRoverRLIons",'d',false);
145 dRoverRPrm2->SetParameterRange("dRoverRLIons>0. && dRoverRLIons<=1.");
146 StepFuncCmd2->SetParameter(dRoverRPrm2);
147
148 G4UIparameter* finalRangePrm2 = new G4UIparameter("finalRangeLIons",'d',false);
149 finalRangePrm2->SetParameterRange("finalRangeLIons>0.");
150 StepFuncCmd2->SetParameter(finalRangePrm2);
151
152 G4UIparameter* unitPrm2 = new G4UIparameter("unit",'s',true);
153 unitPrm2->SetDefaultValue("mm");
154 StepFuncCmd2->SetParameter(unitPrm2);
155
156 StepFuncCmd3 = new G4UIcommand("/process/eLoss/StepFunctionIons",this);
157 StepFuncCmd3->SetGuidance("Set the energy loss step limitation parameters for ions.");
158 StepFuncCmd3->SetGuidance(" dRoverR : max Range variation per step");
159 StepFuncCmd3->SetGuidance(" finalRange: range for final step");
162
163 G4UIparameter* dRoverRPrm3 = new G4UIparameter("dRoverRMuHad",'d',false);
164 dRoverRPrm3->SetParameterRange("dRoverRIons>0. && dRoverRIons<=1.");
165 StepFuncCmd3->SetParameter(dRoverRPrm3);
166
167 G4UIparameter* finalRangePrm3 = new G4UIparameter("finalRangeIons",'d',false);
168 finalRangePrm3->SetParameterRange("finalRangeIons>0.");
169 StepFuncCmd3->SetParameter(finalRangePrm3);
170
171 G4UIparameter* unitPrm3 = new G4UIparameter("unit",'s',true);
172 unitPrm3->SetDefaultValue("mm");
173 StepFuncCmd3->SetParameter(unitPrm3);
174
175 bfCmd = new G4UIcommand("/process/em/setBiasingFactor",this);
176 bfCmd->SetGuidance("Set factor for the process cross section.");
177 bfCmd->SetGuidance(" procName : process name");
178 bfCmd->SetGuidance(" procFact : factor");
179 bfCmd->SetGuidance(" flagFact : flag to change weight");
182
183 G4UIparameter* procName = new G4UIparameter("procName",'s',false);
184 bfCmd->SetParameter(procName);
185
186 G4UIparameter* procFact = new G4UIparameter("procFact",'d',false);
187 bfCmd->SetParameter(procFact);
188
189 G4UIparameter* flagFact = new G4UIparameter("flagFact",'s',false);
190 bfCmd->SetParameter(flagFact);
191
192 fiCmd = new G4UIcommand("/process/em/setForcedInteraction",this);
193 fiCmd->SetGuidance("Set factor for the process cross section.");
194 fiCmd->SetGuidance(" procNam : process name");
195 fiCmd->SetGuidance(" regNam : region name");
196 fiCmd->SetGuidance(" tlength : fixed target length");
197 fiCmd->SetGuidance(" unitT : length unit");
198 fiCmd->SetGuidance(" tflag : flag to change weight");
201
202 G4UIparameter* procNam = new G4UIparameter("procNam",'s',false);
203 fiCmd->SetParameter(procNam);
204
205 G4UIparameter* regNam = new G4UIparameter("regNam",'s',false);
206 fiCmd->SetParameter(regNam);
207
208 G4UIparameter* tlength = new G4UIparameter("tlength",'d',false);
209 tlength->SetParameterRange("tlength>0");
210 fiCmd->SetParameter(tlength);
211
212 G4UIparameter* unitT = new G4UIparameter("unitT",'s',true);
213 unitT->SetDefaultUnit("mm");
214 fiCmd->SetParameter(unitT);
215
216 G4UIparameter* flagT = new G4UIparameter("tflag",'b',true);
217 flagT->SetDefaultValue(true);
218 fiCmd->SetParameter(flagT);
219
220 bsCmd = new G4UIcommand("/process/em/setSecBiasing",this);
221 bsCmd->SetGuidance("Set bremsstrahlung or delta-e- splitting/Russian roulette per region.");
222 bsCmd->SetGuidance(" bProcNam : process name");
223 bsCmd->SetGuidance(" bRegNam : region name");
224 bsCmd->SetGuidance(" bFactor : number of split gamma or probability of Russian roulette");
225 bsCmd->SetGuidance(" bEnergy : max energy of a secondary for this biasing method");
226 bsCmd->SetGuidance(" bUnit : energy unit");
229
230 G4UIparameter* bProcNam = new G4UIparameter("bProcNam",'s',false);
231 bsCmd->SetParameter(bProcNam);
232
233 G4UIparameter* bRegNam = new G4UIparameter("bRegNam",'s',false);
234 bsCmd->SetParameter(bRegNam);
235
236 G4UIparameter* bFactor = new G4UIparameter("bFactor",'d',false);
237 bsCmd->SetParameter(bFactor);
238
239 G4UIparameter* bEnergy = new G4UIparameter("bEnergy",'d',false);
240 bsCmd->SetParameter(bEnergy);
241
242 G4UIparameter* bUnit = new G4UIparameter("bUnit",'s',true);
243 bUnit->SetDefaultUnit("MeV");
244 bsCmd->SetParameter(bUnit);
245
246 dirSplitCmd = new G4UIcmdWithABool("/process/em/setDirectionalSplitting",this);
247 dirSplitCmd->SetGuidance("Enable directional brem splitting");
250
251 qeCmd = new G4UIcmdWithABool("/process/em/QuantumEntanglement",this);
252 qeCmd->SetGuidance("Enable quantum entanglement");
255
256 dirSplitTargetCmd = new G4UIcmdWith3VectorAndUnit("/process/em/setDirectionalSplittingTarget",this);
257 dirSplitTargetCmd->SetGuidance("Position of arget for directional splitting");
259
260 dirSplitRadiusCmd = new G4UIcmdWithADoubleAndUnit("/process/em/setDirectionalSplittingRadius",this);
261 dirSplitRadiusCmd->SetGuidance("Radius of target for directional splitting");
264}
@ G4State_Idle
@ G4State_PreInit
G4UIcmdWith3VectorAndUnit * dirSplitTargetCmd
G4UIcmdWithADoubleAndUnit * dirSplitRadiusCmd
void SetToBeBroadcasted(G4bool val)
Definition: G4UIcommand.hh:172
void SetParameter(G4UIparameter *const newParameter)
Definition: G4UIcommand.hh:146
void SetGuidance(const char *aGuidance)
Definition: G4UIcommand.hh:156
void AvailableForStates(G4ApplicationState s1)
Definition: G4UIcommand.cc:288
void SetDefaultValue(const char *theDefaultValue)
void SetParameterRange(const char *theRange)
void SetParameterCandidates(const char *theString)
void SetDefaultUnit(const char *theDefaultUnit)

References G4UIcommand::AvailableForStates(), bfCmd, bsCmd, dirSplitCmd, dirSplitRadiusCmd, dirSplitTargetCmd, fiCmd, G4State_Idle, G4State_PreInit, mscoCmd, paiCmd, qeCmd, G4UIparameter::SetDefaultUnit(), G4UIparameter::SetDefaultValue(), G4UIcommand::SetGuidance(), G4UIcommand::SetParameter(), G4UIparameter::SetParameterCandidates(), G4UIparameter::SetParameterRange(), G4UIcommand::SetToBeBroadcasted(), StepFuncCmd, StepFuncCmd1, StepFuncCmd2, StepFuncCmd3, and SubSecCmd.

◆ ~G4EmExtraParametersMessenger()

G4EmExtraParametersMessenger::~G4EmExtraParametersMessenger ( )
override

Definition at line 268 of file G4EmExtraParametersMessenger.cc.

269{
270 delete paiCmd;
271 delete mscoCmd;
272 delete SubSecCmd;
273 delete bfCmd;
274 delete fiCmd;
275 delete bsCmd;
276 delete qeCmd;
277 delete StepFuncCmd;
278 delete StepFuncCmd1;
279 delete StepFuncCmd2;
280 delete StepFuncCmd3;
281 delete dirSplitCmd;
282 delete dirSplitTargetCmd;
283 delete dirSplitRadiusCmd;
284}

References bfCmd, bsCmd, dirSplitCmd, dirSplitRadiusCmd, dirSplitTargetCmd, fiCmd, mscoCmd, paiCmd, qeCmd, StepFuncCmd, StepFuncCmd1, StepFuncCmd2, StepFuncCmd3, and SubSecCmd.

◆ G4EmExtraParametersMessenger() [2/2]

G4EmExtraParametersMessenger::G4EmExtraParametersMessenger ( const G4EmExtraParametersMessenger )
delete

Member Function Documentation

◆ AddUIcommand()

void G4UImessenger::AddUIcommand ( G4UIcommand newCommand)
protectedinherited

Definition at line 149 of file G4UImessenger.cc.

150{
151 G4cerr << "Warning : Old style definition of G4UIcommand <"
152 << newCommand->GetCommandPath() << ">." << G4endl;
153}
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition: G4ios.hh:57
const G4String & GetCommandPath() const
Definition: G4UIcommand.hh:136

References G4cerr, G4endl, and G4UIcommand::GetCommandPath().

◆ BtoS()

G4String G4UImessenger::BtoS ( G4bool  b)
protectedinherited

Definition at line 98 of file G4UImessenger.cc.

99{
100 G4String vl = "0";
101 if(b)
102 vl = "true";
103 return vl;
104}

◆ CommandsShouldBeInMaster()

G4bool G4UImessenger::CommandsShouldBeInMaster ( ) const
inlineinherited

Definition at line 77 of file G4UImessenger.hh.

78 {
80 }
G4bool commandsShouldBeInMaster

References G4UImessenger::commandsShouldBeInMaster.

Referenced by G4UIcommand::G4UIcommandCommonConstructorCode().

◆ CreateCommand()

template<typename T >
T * G4UImessenger::CreateCommand ( const G4String cname,
const G4String dsc 
)
protectedinherited

Definition at line 110 of file G4UImessenger.hh.

111{
112 G4String path;
113 if(cname[0] != '/')
114 {
115 path = baseDirName + cname;
116 if(path[0] != '/')
117 path = "/" + path;
118 }
119
120 T* command = new T(path.c_str(), this);
121 command->SetGuidance(dsc.c_str());
122
123 return command;
124}
G4String baseDirName

References G4UImessenger::baseDirName.

◆ CreateDirectory()

void G4UImessenger::CreateDirectory ( const G4String path,
const G4String dsc,
G4bool  commandsToBeBroadcasted = true 
)
protectedinherited

Definition at line 156 of file G4UImessenger.cc.

158{
160
161 G4String fullpath = path;
162 if(fullpath.back() != '/')
163 fullpath.append("/");
164
165 G4UIcommandTree* tree = ui->GetTree()->FindCommandTree(fullpath.c_str());
166 if(tree != nullptr)
167 {
168 baseDirName = tree->GetPathName();
169 }
170 else
171 {
172 baseDir = new G4UIdirectory(fullpath.c_str(), commandsToBeBroadcasted);
173 baseDirName = fullpath;
174 baseDir->SetGuidance(dsc.c_str());
175 }
176}
const G4String & GetPathName() const
G4UIcommandTree * FindCommandTree(const char *commandPath)
G4UIcommandTree * GetTree() const
Definition: G4UImanager.hh:186
static G4UImanager * GetUIpointer()
Definition: G4UImanager.cc:77
G4UIdirectory * baseDir

References G4UImessenger::baseDir, G4UImessenger::baseDirName, G4UIcommandTree::FindCommandTree(), G4UIcommandTree::GetPathName(), G4UImanager::GetTree(), G4UImanager::GetUIpointer(), and G4UIcommand::SetGuidance().

Referenced by G4MoleculeShootMessenger::G4MoleculeShootMessenger(), and G4UImessenger::G4UImessenger().

◆ DtoS()

G4String G4UImessenger::DtoS ( G4double  a)
protectedinherited

Definition at line 90 of file G4UImessenger.cc.

91{
92 std::ostringstream os;
93 os << a;
94 return G4String(os.str());
95}

Referenced by G4ScoreQuantityMessenger::FilterCommands(), and G4UIcontrolMessenger::SetNewValue().

◆ GetCurrentValue()

G4String G4UImessenger::GetCurrentValue ( G4UIcommand command)
virtualinherited

Reimplemented in G4ScoreQuantityMessenger, G4VisCommandModelCreate< Factory >, G4VisCommandListManagerList< Manager >, G4VisCommandListManagerSelect< Manager >, G4VisCommandManagerMode< Manager >, G4ToolsAnalysisMessenger, G4ScoringMessenger, G4EvManMessenger, G4GeneralParticleSourceMessenger, G4ParticleGunMessenger, G4GeometryMessenger, G4GenericMessenger, G4UIcontrolMessenger, GFlashShowerModelMessenger, G4DecayTableMessenger, G4ParticleMessenger, G4ParticlePropertyMessenger, G4tgrMessenger, G4PersistencyCenterMessenger, G4ProductionCutsTableMessenger, G4SchedulerMessenger, G4VITSteppingVerbose, G4MoleculeShootMessenger, G4MoleculeGunMessenger, G4ProcessManagerMessenger, G4ProcessTableMessenger, G4MatScanMessenger, G4RunMessenger, G4UserPhysicsListMessenger, G4TrackingMessenger, G4GMocrenMessenger, G4HepRepMessenger, G4VisCommandAbortReviewKeptEvents, G4VisCommandDrawOnlyToBeKeptEvents, G4VisCommandEnable, G4VisCommandList, G4VisCommandReviewKeptEvents, G4VisCommandVerbose, G4VisCommandGeometryList, G4VisCommandGeometryRestore, G4VisCommandGeometrySetColour, G4VisCommandGeometrySetDaughtersInvisible, G4VisCommandGeometrySetForceAuxEdgeVisible, G4VisCommandGeometrySetForceCloud, G4VisCommandGeometrySetForceSolid, G4VisCommandGeometrySetForceLineSegmentsPerCircle, G4VisCommandGeometrySetForceWireframe, G4VisCommandGeometrySetLineStyle, G4VisCommandGeometrySetLineWidth, G4VisCommandGeometrySetVisibility, G4VisCommandSceneActivateModel, G4VisCommandSceneCreate, G4VisCommandSceneEndOfEventAction, G4VisCommandSceneEndOfRunAction, G4VisCommandSceneList, G4VisCommandSceneNotifyHandlers, G4VisCommandSceneRemoveModel, G4VisCommandSceneSelect, G4VisCommandSceneShowExtents, G4VisCommandSceneAddArrow, G4VisCommandSceneAddArrow2D, G4VisCommandSceneAddAxes, G4VisCommandSceneAddDate, G4VisCommandSceneAddDigis, G4VisCommandSceneAddEventID, G4VisCommandSceneAddExtent, G4VisCommandSceneAddElectricField, G4VisCommandSceneAddFrame, G4VisCommandSceneAddGPS, G4VisCommandSceneAddGhosts, G4VisCommandSceneAddHits, G4VisCommandSceneAddLine, G4VisCommandSceneAddLine2D, G4VisCommandSceneAddLocalAxes, G4VisCommandSceneAddLogicalVolume, G4VisCommandSceneAddLogo, G4VisCommandSceneAddLogo2D, G4VisCommandSceneAddMagneticField, G4VisCommandSceneAddPSHits, G4VisCommandSceneAddScale, G4VisCommandSceneAddText, G4VisCommandSceneAddText2D, G4VisCommandSceneAddTrajectories, G4VisCommandSceneAddUserAction, G4VisCommandSceneAddVolume, G4VisCommandSceneAddPlotter, G4VisCommandSceneHandlerAttach, G4VisCommandSceneHandlerCreate, G4VisCommandSceneHandlerList, G4VisCommandSceneHandlerSelect, G4VisCommandSetArrow3DLineSegmentsPerCircle, G4VisCommandSetColour, G4VisCommandSetExtentForField, G4VisCommandSetLineWidth, G4VisCommandSetTextColour, G4VisCommandSetTextLayout, G4VisCommandSetTextSize, G4VisCommandSetTouchable, G4VisCommandSetVolumeForField, G4VisCommandsTouchable, G4VisCommandsTouchableSet, G4VisCommandViewerAddCutawayPlane, G4VisCommandViewerCentreOn, G4VisCommandViewerChangeCutawayPlane, G4VisCommandViewerClear, G4VisCommandViewerClearCutawayPlanes, G4VisCommandViewerClearTransients, G4VisCommandViewerClearVisAttributesModifiers, G4VisCommandViewerClone, G4VisCommandViewerColourByDensity, G4VisCommandViewerCopyViewFrom, G4VisCommandViewerCreate, G4VisCommandViewerDolly, G4VisCommandViewerFlush, G4VisCommandViewerInterpolate, G4VisCommandViewerList, G4VisCommandViewerPan, G4VisCommandViewerReset, G4VisCommandViewerRefresh, G4VisCommandViewerRebuild, G4VisCommandViewerSave, G4VisCommandViewerScale, G4VisCommandViewerSelect, G4VisCommandViewerUpdate, G4VisCommandViewerZoom, G4VisCommandViewerDefaultHiddenEdge, G4VisCommandViewerDefaultStyle, G4VisCommandsViewerSet, G4VModelCommand< T >, G4VModelCommand< M >, G4RTMessenger, G4ASCIITreeMessenger, G4VtkMessenger, G4PolarizationMessenger, and G4DNAChemistryManager.

Definition at line 58 of file G4UImessenger.cc.

59{
60 G4String nullString;
61 return nullString;
62}

Referenced by G4UIcommand::DoIt(), and G4UIcommand::GetCurrentValue().

◆ ItoS()

G4String G4UImessenger::ItoS ( G4int  i)
protectedinherited

Definition at line 82 of file G4UImessenger.cc.

83{
84 std::ostringstream os;
85 os << i;
86 return G4String(os.str());
87}

Referenced by G4GenericMessenger::DeclareMethod(), and G4ParticleGunMessenger::GetCurrentValue().

◆ operator!=()

G4bool G4UImessenger::operator!= ( const G4UImessenger messenger) const
inherited

Definition at line 76 of file G4UImessenger.cc.

77{
78 return this != &messenger;
79}

◆ operator=()

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

◆ operator==()

G4bool G4UImessenger::operator== ( const G4UImessenger messenger) const
inherited

Definition at line 70 of file G4UImessenger.cc.

71{
72 return this == &messenger;
73}

◆ SetNewValue()

void G4EmExtraParametersMessenger::SetNewValue ( G4UIcommand command,
G4String  newValue 
)
overridevirtual

Reimplemented from G4UImessenger.

Definition at line 288 of file G4EmExtraParametersMessenger.cc.

290{
291 G4bool physicsModified = false;
292
293 if (command == paiCmd) {
294 G4String s1(""),s2(""),s3("");
295 std::istringstream is(newValue);
296 is >> s1 >> s2 >> s3;
297 theParameters->AddPAIModel(s1, s2, s3);
298 } else if (command == mscoCmd) {
299 G4String s1(""),s2("");
300 std::istringstream is(newValue);
301 is >> s1 >> s2;
302 theParameters->AddPhysics(s1, s2);
303 } else if (command == StepFuncCmd || command == StepFuncCmd1 || command == StepFuncCmd2 || command == StepFuncCmd3) {
304 G4double v1,v2;
305 G4String unt;
306 std::istringstream is(newValue);
307 is >> v1 >> v2 >> unt;
308 v2 *= G4UIcommand::ValueOf(unt);
309 if(command == StepFuncCmd) {
311 } else if(command == StepFuncCmd1) {
313 } else if(command == StepFuncCmd2) {
315 } else {
317 }
318 physicsModified = true;
319 } else if (command == SubSecCmd) {
321 } else if (command == bfCmd) {
322 G4double v1(1.0);
323 G4String s0(""),s1("");
324 std::istringstream is(newValue);
325 is >> s0 >> v1 >> s1;
326 G4bool yes = false;
327 if(s1 == "true") { yes = true; }
329 physicsModified = true;
330 } else if (command == fiCmd) {
331 G4double v1(0.0);
332 G4String s1(""),s2(""),s3(""),unt("mm");
333 std::istringstream is(newValue);
334 is >> s1 >> s2 >> v1 >> unt >> s3;
335 G4bool yes = false;
336 if(s3 == "true") { yes = true; }
337 v1 *= G4UIcommand::ValueOf(unt);
339 physicsModified = true;
340 } else if (command == bsCmd) {
341 G4double fb(1.0),en(1.e+30);
342 G4String s1(""),s2(""),unt("MeV");
343 std::istringstream is(newValue);
344 is >> s1 >> s2 >> fb >> en >> unt;
345 en *= G4UIcommand::ValueOf(unt);
347 physicsModified = true;
348 } else if (command == qeCmd) {
350 } else if (command == dirSplitCmd) {
352 dirSplitCmd->GetNewBoolValue(newValue));
353 physicsModified = true;
354 } else if (command == dirSplitTargetCmd) {
357 physicsModified = true;
358 } else if (command == dirSplitRadiusCmd) {
361 physicsModified = true;
362 }
363
364 if(physicsModified) {
365 G4UImanager::GetUIpointer()->ApplyCommand("/run/physicsModified");
366 }
367}
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
void SetQuantumEntanglement(G4bool v)
void SetStepFunctionMuHad(G4double v1, G4double v2)
void SetDirectionalSplittingRadius(G4double r)
void SetStepFunctionIons(G4double v1, G4double v2)
void ActivateForcedInteraction(const G4String &procname, const G4String &region, G4double length, G4bool wflag)
void ActivateSecondaryBiasing(const G4String &name, const G4String &region, G4double factor, G4double energyLimit)
void SetDirectionalSplitting(G4bool v)
void SetDirectionalSplittingTarget(const G4ThreeVector &v)
void SetSubCutRegion(const G4String &region)
void SetProcessBiasingFactor(const G4String &procname, G4double val, G4bool wflag)
void SetStepFunction(G4double v1, G4double v2)
void SetStepFunctionLightIons(G4double v1, G4double v2)
void AddPAIModel(const G4String &particle, const G4String &region, const G4String &type)
void AddPhysics(const G4String &region, const G4String &type)
static G4ThreeVector GetNew3VectorValue(const char *paramString)
static G4bool GetNewBoolValue(const char *paramString)
static G4double GetNewDoubleValue(const char *paramString)
static G4double ValueOf(const char *unitName)
Definition: G4UIcommand.cc:363
G4int ApplyCommand(const char *aCommand)
Definition: G4UImanager.cc:485

References G4EmExtraParameters::ActivateForcedInteraction(), G4EmExtraParameters::ActivateSecondaryBiasing(), G4EmExtraParameters::AddPAIModel(), G4EmExtraParameters::AddPhysics(), G4UImanager::ApplyCommand(), bfCmd, bsCmd, dirSplitCmd, dirSplitRadiusCmd, dirSplitTargetCmd, fiCmd, G4UIcmdWith3VectorAndUnit::GetNew3VectorValue(), G4UIcmdWithABool::GetNewBoolValue(), G4UIcmdWithADoubleAndUnit::GetNewDoubleValue(), G4UImanager::GetUIpointer(), mscoCmd, paiCmd, qeCmd, G4InuclParticleNames::s0, G4EmExtraParameters::SetDirectionalSplitting(), G4EmExtraParameters::SetDirectionalSplittingRadius(), G4EmExtraParameters::SetDirectionalSplittingTarget(), G4EmExtraParameters::SetProcessBiasingFactor(), G4EmExtraParameters::SetQuantumEntanglement(), G4EmExtraParameters::SetStepFunction(), G4EmExtraParameters::SetStepFunctionIons(), G4EmExtraParameters::SetStepFunctionLightIons(), G4EmExtraParameters::SetStepFunctionMuHad(), G4EmExtraParameters::SetSubCutRegion(), StepFuncCmd, StepFuncCmd1, StepFuncCmd2, StepFuncCmd3, SubSecCmd, theParameters, and G4UIcommand::ValueOf().

◆ StoB()

G4bool G4UImessenger::StoB ( G4String  s)
protectedinherited

Definition at line 137 of file G4UImessenger.cc.

138{
140 G4bool vl = false;
141 if(v == "Y" || v == "YES" || v == "1" || v == "T" || v == "TRUE")
142 {
143 vl = true;
144 }
145 return vl;
146}
G4String to_upper_copy(G4String str)
Return uppercase copy of string.

References G4StrUtil::to_upper_copy().

Referenced by G4LocalThreadCoutMessenger::SetNewValue(), G4CascadeParamMessenger::SetNewValue(), G4ScoreQuantityMessenger::SetNewValue(), and G4ScoringMessenger::SetNewValue().

◆ StoD()

G4double G4UImessenger::StoD ( G4String  s)
protectedinherited

◆ StoI()

G4int G4UImessenger::StoI ( G4String  s)
protectedinherited

◆ StoL()

G4long G4UImessenger::StoL ( G4String  s)
protectedinherited

Definition at line 117 of file G4UImessenger.cc.

118{
119 G4long vl;
120 const char* t = str;
121 std::istringstream is(t);
122 is >> vl;
123 return vl;
124}
long G4long
Definition: G4Types.hh:87

Referenced by G4RunMessenger::SetNewValue().

Field Documentation

◆ baseDir

G4UIdirectory* G4UImessenger::baseDir = nullptr
protectedinherited

◆ baseDirName

G4String G4UImessenger::baseDirName = ""
protectedinherited

◆ bfCmd

G4UIcommand* G4EmExtraParametersMessenger::bfCmd
private

◆ bsCmd

G4UIcommand* G4EmExtraParametersMessenger::bsCmd
private

◆ commandsShouldBeInMaster

G4bool G4UImessenger::commandsShouldBeInMaster = false
protectedinherited

◆ dirSplitCmd

G4UIcmdWithABool* G4EmExtraParametersMessenger::dirSplitCmd
private

◆ dirSplitRadiusCmd

G4UIcmdWithADoubleAndUnit* G4EmExtraParametersMessenger::dirSplitRadiusCmd
private

◆ dirSplitTargetCmd

G4UIcmdWith3VectorAndUnit* G4EmExtraParametersMessenger::dirSplitTargetCmd
private

◆ fiCmd

G4UIcommand* G4EmExtraParametersMessenger::fiCmd
private

◆ mscoCmd

G4UIcommand* G4EmExtraParametersMessenger::mscoCmd
private

◆ paiCmd

G4UIcommand* G4EmExtraParametersMessenger::paiCmd
private

◆ qeCmd

G4UIcmdWithABool* G4EmExtraParametersMessenger::qeCmd
private

◆ StepFuncCmd

G4UIcommand* G4EmExtraParametersMessenger::StepFuncCmd
private

◆ StepFuncCmd1

G4UIcommand* G4EmExtraParametersMessenger::StepFuncCmd1
private

◆ StepFuncCmd2

G4UIcommand* G4EmExtraParametersMessenger::StepFuncCmd2
private

◆ StepFuncCmd3

G4UIcommand* G4EmExtraParametersMessenger::StepFuncCmd3
private

◆ SubSecCmd

G4UIcmdWithAString* G4EmExtraParametersMessenger::SubSecCmd
private

◆ theParameters

G4EmExtraParameters* G4EmExtraParametersMessenger::theParameters
private

Definition at line 79 of file G4EmExtraParametersMessenger.hh.

Referenced by SetNewValue().


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