Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4FieldParametersMessenger.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// Implementation of the G4FieldParametersMessenger class
27//
28// Author: Ivana Hrivnacova (IJClab, Orsay), 2024.
29// --------------------------------------------------------------------
30
32#include "G4FieldParameters.hh"
33
34#include "G4UIcmdWithABool.hh"
35#include "G4UIcmdWithADouble.hh"
37#include "G4UIcmdWithAString.hh"
39#include "G4UIdirectory.hh"
40
41//_____________________________________________________________________________
43 G4FieldParameters* fieldParameters)
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}
191
192//_____________________________________________________________________________
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}
209
210//
211// public methods
212//
213
214//_____________________________________________________________________________
216 G4UIcommand* command, G4String newValues)
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}
@ 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
void SetNewValue(G4UIcommand *command, G4String newValues) override
G4FieldParameters defines the type of equation of motion of a particle in a field and the integration...
const G4String & GetVolumeName() const
static G4String EquationTypeName(G4EquationType equation)
static G4String StepperTypeName(G4StepperType stepper)
static G4String FieldTypeName(G4FieldType field)