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

G4FieldParametersMessenger is a messenger class that defines commands for field configuration. Only equation type and stepper type values that are handled by G4FieldBuilder are accepted by the commands. More...

#include <G4FieldParametersMessenger.hh>

Inheritance diagram for G4FieldParametersMessenger:

Public Member Functions

 G4FieldParametersMessenger (G4FieldParameters *fieldParameters)
 ~G4FieldParametersMessenger () override
 G4FieldParametersMessenger ()=delete
 G4FieldParametersMessenger (const G4FieldParametersMessenger &)=delete
G4FieldParametersMessengeroperator= (const G4FieldParametersMessenger &)=delete
void SetNewValue (G4UIcommand *command, G4String newValues) override
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

G4FieldParametersMessenger is a messenger class that defines commands for field configuration. Only equation type and stepper type values that are handled by G4FieldBuilder are accepted by the commands.

Definition at line 80 of file G4FieldParametersMessenger.hh.

Constructor & Destructor Documentation

◆ G4FieldParametersMessenger() [1/3]

G4FieldParametersMessenger::G4FieldParametersMessenger ( G4FieldParameters * fieldParameters)

Standard constructor for G4FieldParametersMessenger.

Parameters
[in]fieldParametersPointer to the field parameters object.

Definition at line 42 of file G4FieldParametersMessenger.cc.

44 : fFieldParameters(fieldParameters)
45{
46 // Standard constructor
47
48 G4String directoryName = "/field/";
49 if (!fieldParameters->GetVolumeName().empty())
50 {
51 directoryName.append(fieldParameters->GetVolumeName());
52 directoryName.append("/");
53 fDirectory = new G4UIdirectory(directoryName);
54 fDirectory->SetGuidance("Magnetic field control commands.");
55 }
56
57 G4String commandName = directoryName;
58 commandName.append("fieldType");
59 fFieldTypeCmd = new G4UIcmdWithAString(commandName, this);
60 G4String guidance = "Select type of the field";
61 fFieldTypeCmd->SetGuidance(guidance);
62 fFieldTypeCmd->SetParameterName("FieldType", false);
63 G4String candidates;
64 for (G4int i = kMagnetic; i <= kGravity; ++i)
65 {
66 auto ft = (G4FieldType)i;
67 candidates += G4FieldParameters::FieldTypeName(ft);
68 candidates += " ";
69 }
70 fFieldTypeCmd->SetCandidates(candidates);
71 fFieldTypeCmd->AvailableForStates(
73
74 commandName = directoryName;
75 commandName.append("equationType");
76 fEquationTypeCmd = new G4UIcmdWithAString(commandName, this);
77 guidance = "Select type of the equation of motion of a particle in a field";
78 fEquationTypeCmd->SetGuidance(guidance);
79 fEquationTypeCmd->SetParameterName("EquationType", false);
80 candidates = "";
81 for (G4int i = kEqMagnetic; i <= kEqEMfieldWithEDM; ++i)
82 {
83 auto et = (G4EquationType)i;
85 candidates += " ";
86 }
87 fEquationTypeCmd->SetCandidates(candidates);
88 fEquationTypeCmd->AvailableForStates(
90
91 commandName = directoryName;
92 commandName.append("stepperType");
93 fStepperTypeCmd = new G4UIcmdWithAString(commandName, this);
94 guidance =
95 "Select type of the the integrator of particle's equation of motion in a "
96 "field";
97 fStepperTypeCmd->SetGuidance(guidance);
98 fStepperTypeCmd->SetParameterName("StepperType", false);
99 candidates = "";
100 for (G4int i = kCashKarpRKF45; i <= kRK547FEq3; ++i)
101 {
102 auto st = (G4StepperType)i;
103 if (st == kUserStepper) { continue; }
104 candidates += G4FieldParameters::StepperTypeName(st);
105 candidates += " ";
106 }
107 fStepperTypeCmd->SetCandidates(candidates);
108 fStepperTypeCmd->AvailableForStates(
110
111 commandName = directoryName;
112 commandName.append("setMinimumStep");
113 fSetMinimumStepCmd = new G4UIcmdWithADoubleAndUnit(commandName, this);
114 fSetMinimumStepCmd->SetGuidance("Set minimum step in G4ChordFinder");
115 fSetMinimumStepCmd->SetParameterName("StepMinimum", false);
116 fSetMinimumStepCmd->SetDefaultUnit("mm");
117 fSetMinimumStepCmd->SetUnitCategory("Length");
118 fSetMinimumStepCmd->AvailableForStates(
120
121 commandName = directoryName;
122 commandName.append("setDeltaChord");
123 fSetDeltaChordCmd = new G4UIcmdWithADoubleAndUnit(commandName, this);
124 fSetDeltaChordCmd->SetGuidance("Set delta chord in G4ChordFinder");
125 fSetDeltaChordCmd->SetParameterName("DeltaChord", false);
126 fSetDeltaChordCmd->SetDefaultUnit("mm");
127 fSetDeltaChordCmd->SetUnitCategory("Length");
128 fSetDeltaChordCmd->AvailableForStates(
130
131 commandName = directoryName;
132 commandName.append("setDeltaOneStep");
133 fSetDeltaOneStepCmd = new G4UIcmdWithADoubleAndUnit(commandName, this);
134 fSetDeltaOneStepCmd->SetGuidance(
135 "Set delta one step in global field manager");
136 fSetDeltaOneStepCmd->SetParameterName("DeltaOneStep", false);
137 fSetDeltaOneStepCmd->SetDefaultUnit("mm");
138 fSetDeltaOneStepCmd->SetUnitCategory("Length");
139 fSetDeltaOneStepCmd->AvailableForStates(
141
142 commandName = directoryName;
143 commandName.append("setDeltaIntersection");
144 fSetDeltaIntersectionCmd = new G4UIcmdWithADoubleAndUnit(commandName, this);
145 fSetDeltaIntersectionCmd->SetGuidance(
146 "Set delta intersection in global field manager");
147 fSetDeltaIntersectionCmd->SetParameterName("DeltaIntersection", false);
148 fSetDeltaIntersectionCmd->SetDefaultUnit("mm");
149 fSetDeltaIntersectionCmd->SetUnitCategory("Length");
150 fSetDeltaIntersectionCmd->AvailableForStates(
152
153 commandName = directoryName;
154 commandName.append("setMinimumEpsilonStep");
155 fSetMinimumEpsilonStepCmd = new G4UIcmdWithADouble(commandName, this);
156 fSetMinimumEpsilonStepCmd->SetGuidance(
157 "Set minimum epsilon step in global field manager");
158 fSetMinimumEpsilonStepCmd->SetParameterName("MinimumEpsilonStep", false);
159 fSetMinimumEpsilonStepCmd->AvailableForStates(
161
162 commandName = directoryName;
163 commandName.append("setMaximumEpsilonStep");
164 fSetMaximumEpsilonStepCmd = new G4UIcmdWithADouble(commandName, this);
165 fSetMaximumEpsilonStepCmd->SetGuidance(
166 "Set maximum epsilon step in global field manager");
167 fSetMaximumEpsilonStepCmd->SetParameterName("MaximumEpsilonStep", false);
168 fSetMaximumEpsilonStepCmd->AvailableForStates(
170
171 commandName = directoryName;
172 commandName.append("setConstDistance");
173 fSetConstDistanceCmd = new G4UIcmdWithADoubleAndUnit(commandName, this);
174 fSetConstDistanceCmd->SetGuidance(
175 "Set the distance within which the field is considered constant.");
176 fSetConstDistanceCmd->SetGuidance(
177 "Non-zero value will trigger creating a cached magnetic field.");
178 fSetConstDistanceCmd->SetParameterName("ConstDistance", false);
179 fSetConstDistanceCmd->SetDefaultUnit("mm");
180 fSetConstDistanceCmd->SetUnitCategory("Length");
181 fSetConstDistanceCmd->SetRange("ConstDistance >= 0");
182 fSetConstDistanceCmd->AvailableForStates(G4State_PreInit);
183
184 commandName = std::move(directoryName);
185 commandName.append("printParameters");
186 fPrintParametersCmd = new G4UIcmdWithoutParameter(commandName, this);
187 fPrintParametersCmd->SetGuidance("Prints all accuracy parameters.");
188 fPrintParametersCmd->AvailableForStates(
190}
@ G4State_Init
@ G4State_Idle
@ G4State_PreInit
G4EquationType
G4EquationType defines the types of equations of motion of a particle in a field in Geant4.
@ kEqMagnetic
@ kEqEMfieldWithEDM
G4FieldType
G4FieldType defines the available fields in Geant4.
@ kGravity
gravity field
@ kMagnetic
magnetic field
G4StepperType
G4StepperType defines the available integrator of particle's equation of motion in Geant4.
@ kCashKarpRKF45
G4CashKarpRKF45.
@ kUserStepper
User defined stepper.
@ kRK547FEq3
G4RK547FEq3.
int G4int
Definition G4Types.hh:85
const G4String & GetVolumeName() const
static G4String EquationTypeName(G4EquationType equation)
static G4String StepperTypeName(G4StepperType stepper)
static G4String FieldTypeName(G4FieldType field)

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

◆ ~G4FieldParametersMessenger()

G4FieldParametersMessenger::~G4FieldParametersMessenger ( )
override

Destructor.

Definition at line 193 of file G4FieldParametersMessenger.cc.

194{
195 // Destructor
196
197 delete fDirectory;
198 delete fFieldTypeCmd;
199 delete fEquationTypeCmd;
200 delete fStepperTypeCmd;
201 delete fSetMinimumStepCmd;
202 delete fSetDeltaChordCmd;
203 delete fSetDeltaOneStepCmd;
204 delete fSetDeltaIntersectionCmd;
205 delete fSetMinimumEpsilonStepCmd;
206 delete fSetMaximumEpsilonStepCmd;
207 delete fSetConstDistanceCmd;
208}

◆ G4FieldParametersMessenger() [2/3]

G4FieldParametersMessenger::G4FieldParametersMessenger ( )
delete

Default constructor, copy constructor and assignment operator not allowed.

◆ G4FieldParametersMessenger() [3/3]

G4FieldParametersMessenger::G4FieldParametersMessenger ( const G4FieldParametersMessenger & )
delete

Member Function Documentation

◆ operator=()

G4FieldParametersMessenger & G4FieldParametersMessenger::operator= ( const G4FieldParametersMessenger & )
delete

◆ SetNewValue()

void G4FieldParametersMessenger::SetNewValue ( G4UIcommand * command,
G4String newValues )
overridevirtual

Applies command to the associated object.

Reimplemented from G4UImessenger.

Definition at line 215 of file G4FieldParametersMessenger.cc.

217{
218 // Apply command to the associated object.
219
220 if (command == fFieldTypeCmd)
221 {
222 for (G4int i = kMagnetic; i <= kGravity; ++i)
223 {
224 auto ft = (G4FieldType)i;
225 if (newValues == G4FieldParameters::FieldTypeName(ft))
226 {
227 fFieldParameters->SetFieldType(ft);
228 break;
229 }
230 }
231 return;
232 }
233
234 if (command == fEquationTypeCmd)
235 {
236 for (G4int i = kEqMagnetic; i <= kEqEMfieldWithEDM; ++i)
237 {
238 auto et = (G4EquationType)i;
239 if (newValues == G4FieldParameters::EquationTypeName(et))
240 {
241 fFieldParameters->SetEquationType(et);
242 break;
243 }
244 }
245 return;
246 }
247
248 if (command == fStepperTypeCmd)
249 {
250 for (G4int i = kCashKarpRKF45; i <= kRK547FEq3; ++i)
251 {
252 auto st = (G4StepperType)i;
253 if (newValues == G4FieldParameters::StepperTypeName(st))
254 {
255 fFieldParameters->SetStepperType(st);
256 break;
257 }
258 }
259 return;
260 }
261
262 if (command == fSetMinimumStepCmd)
263 {
264 fFieldParameters->SetMinimumStep(
265 fSetMinimumStepCmd->GetNewDoubleValue(newValues));
266 return;
267 }
268
269 if (command == fSetDeltaChordCmd)
270 {
271 fFieldParameters->SetDeltaChord(
272 fSetDeltaChordCmd->GetNewDoubleValue(newValues));
273 return;
274 }
275
276 if (command == fSetDeltaOneStepCmd)
277 {
278 fFieldParameters->SetDeltaOneStep(
279 fSetDeltaOneStepCmd->GetNewDoubleValue(newValues));
280 return;
281 }
282
283 if (command == fSetDeltaIntersectionCmd)
284 {
285 fFieldParameters->SetDeltaIntersection(
286 fSetDeltaIntersectionCmd->GetNewDoubleValue(newValues));
287 return;
288 }
289
290 if (command == fSetMinimumEpsilonStepCmd)
291 {
292 fFieldParameters->SetMinimumEpsilonStep(
293 fSetMinimumEpsilonStepCmd->GetNewDoubleValue(newValues));
294 return;
295 }
296
297 if (command == fSetMaximumEpsilonStepCmd)
298 {
299 fFieldParameters->SetMaximumEpsilonStep(
300 fSetMaximumEpsilonStepCmd->GetNewDoubleValue(newValues));
301 return;
302 }
303
304 if (command == fSetConstDistanceCmd)
305 {
306 fFieldParameters->SetConstDistance(
307 fSetConstDistanceCmd->GetNewDoubleValue(newValues));
308 return;
309 }
310
311 if (command == fPrintParametersCmd)
312 {
313 fFieldParameters->PrintParameters();
314 return;
315 }
316}

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