Geant4 11.3.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
27/// \file G4FieldParametersMessenger.cc
28/// \brief Implementation of the G4FieldParametersMessenger class
29///
30/// This code was initially developed in Geant4 VMC package
31/// (https://github.com/vmc-project)
32/// and adapted to Geant4.
33///
34/// \author I. Hrivnacova; IJCLab, Orsay
35
37#include "G4FieldParameters.hh"
38
39#include "G4UIcmdWithABool.hh"
40#include "G4UIcmdWithADouble.hh"
42#include "G4UIcmdWithAString.hh"
44#include "G4UIdirectory.hh"
45
46//_____________________________________________________________________________
48 G4FieldParameters* fieldParameters)
49 : fFieldParameters(fieldParameters)
50{
51 // Standard constructor
52
53 G4String directoryName = "/field/";
54 if (fieldParameters->GetVolumeName() != "") {
55 directoryName.append(fieldParameters->GetVolumeName());
56 directoryName.append("/");
57 fDirectory = new G4UIdirectory(directoryName);
58 fDirectory->SetGuidance("Magnetic field control commands.");
59 }
60
61 G4String commandName = directoryName;
62 commandName.append("fieldType");
63 fFieldTypeCmd = new G4UIcmdWithAString(commandName, this);
64 G4String guidance = "Select type of the field";
65 fFieldTypeCmd->SetGuidance(guidance);
66 fFieldTypeCmd->SetParameterName("FieldType", false);
67 G4String candidates;
68 for (G4int i = kMagnetic; i <= kGravity; i++) {
70 candidates += G4FieldParameters::FieldTypeName(ft);
71 candidates += " ";
72 }
73 fFieldTypeCmd->SetCandidates(candidates);
74 fFieldTypeCmd->AvailableForStates(
76
77 commandName = directoryName;
78 commandName.append("equationType");
79 fEquationTypeCmd = new G4UIcmdWithAString(commandName, this);
80 guidance = "Select type of the equation of motion of a particle in a field";
81 fEquationTypeCmd->SetGuidance(guidance);
82 fEquationTypeCmd->SetParameterName("EquationType", false);
83 candidates = "";
84 for (G4int i = kEqMagnetic; i <= kEqEMfieldWithEDM; i++) {
87 candidates += " ";
88 }
89 fEquationTypeCmd->SetCandidates(candidates);
90 fEquationTypeCmd->AvailableForStates(
92
93 commandName = directoryName;
94 commandName.append("stepperType");
95 fStepperTypeCmd = new G4UIcmdWithAString(commandName, this);
96 guidance =
97 "Select type of the the integrator of particle's equation of motion in a "
98 "field";
99 fStepperTypeCmd->SetGuidance(guidance);
100 fStepperTypeCmd->SetParameterName("StepperType", false);
101 candidates = "";
102 for (G4int i = kCashKarpRKF45; i <= kRK547FEq3; i++) {
104 if (st == kUserStepper) continue;
105 candidates += G4FieldParameters::StepperTypeName(st);
106 candidates += " ";
107 }
108 fStepperTypeCmd->SetCandidates(candidates);
109 fStepperTypeCmd->AvailableForStates(
111
112 commandName = directoryName;
113 commandName.append("setMinimumStep");
114 fSetMinimumStepCmd = new G4UIcmdWithADoubleAndUnit(commandName, this);
115 fSetMinimumStepCmd->SetGuidance("Set minimum step in G4ChordFinder");
116 fSetMinimumStepCmd->SetParameterName("StepMinimum", false);
117 fSetMinimumStepCmd->SetDefaultUnit("mm");
118 fSetMinimumStepCmd->SetUnitCategory("Length");
119 fSetMinimumStepCmd->AvailableForStates(
121
122 commandName = directoryName;
123 commandName.append("setDeltaChord");
124 fSetDeltaChordCmd = new G4UIcmdWithADoubleAndUnit(commandName, this);
125 fSetDeltaChordCmd->SetGuidance("Set delta chord in G4ChordFinder");
126 fSetDeltaChordCmd->SetParameterName("DeltaChord", false);
127 fSetDeltaChordCmd->SetDefaultUnit("mm");
128 fSetDeltaChordCmd->SetUnitCategory("Length");
129 fSetDeltaChordCmd->AvailableForStates(
131
132 commandName = directoryName;
133 commandName.append("setDeltaOneStep");
134 fSetDeltaOneStepCmd = new G4UIcmdWithADoubleAndUnit(commandName, this);
135 fSetDeltaOneStepCmd->SetGuidance(
136 "Set delta one step in global field manager");
137 fSetDeltaOneStepCmd->SetParameterName("DeltaOneStep", false);
138 fSetDeltaOneStepCmd->SetDefaultUnit("mm");
139 fSetDeltaOneStepCmd->SetUnitCategory("Length");
140 fSetDeltaOneStepCmd->AvailableForStates(
142
143 commandName = directoryName;
144 commandName.append("setDeltaIntersection");
145 fSetDeltaIntersectionCmd = new G4UIcmdWithADoubleAndUnit(commandName, this);
146 fSetDeltaIntersectionCmd->SetGuidance(
147 "Set delta intersection in global field manager");
148 fSetDeltaIntersectionCmd->SetParameterName("DeltaIntersection", false);
149 fSetDeltaIntersectionCmd->SetDefaultUnit("mm");
150 fSetDeltaIntersectionCmd->SetUnitCategory("Length");
151 fSetDeltaIntersectionCmd->AvailableForStates(
153
154 commandName = directoryName;
155 commandName.append("setMinimumEpsilonStep");
156 fSetMinimumEpsilonStepCmd = new G4UIcmdWithADouble(commandName, this);
157 fSetMinimumEpsilonStepCmd->SetGuidance(
158 "Set minimum epsilon step in global field manager");
159 fSetMinimumEpsilonStepCmd->SetParameterName("MinimumEpsilonStep", false);
160 fSetMinimumEpsilonStepCmd->AvailableForStates(
162
163 commandName = directoryName;
164 commandName.append("setMaximumEpsilonStep");
165 fSetMaximumEpsilonStepCmd = new G4UIcmdWithADouble(commandName, this);
166 fSetMaximumEpsilonStepCmd->SetGuidance(
167 "Set maximum epsilon step in global field manager");
168 fSetMaximumEpsilonStepCmd->SetParameterName("MaximumEpsilonStep", false);
169 fSetMaximumEpsilonStepCmd->AvailableForStates(
171
172 commandName = directoryName;
173 commandName.append("setConstDistance");
174 fSetConstDistanceCmd = new G4UIcmdWithADoubleAndUnit(commandName, this);
175 fSetConstDistanceCmd->SetGuidance(
176 "Set the distance within which the field is considered constant.");
177 fSetConstDistanceCmd->SetGuidance(
178 "Non-zero value will trigger creating a cached magnetic field.");
179 fSetConstDistanceCmd->SetParameterName("ConstDistance", false);
180 fSetConstDistanceCmd->SetDefaultUnit("mm");
181 fSetConstDistanceCmd->SetUnitCategory("Length");
182 fSetConstDistanceCmd->SetRange("ConstDistance >= 0");
183 fSetConstDistanceCmd->AvailableForStates(G4State_PreInit);
184
185 commandName = std::move(directoryName);
186 commandName.append("printParameters");
187 fPrintParametersCmd = new G4UIcmdWithoutParameter(commandName, this);
188 fPrintParametersCmd->SetGuidance("Prints all accuracy parameters.");
189 fPrintParametersCmd->AvailableForStates(
191}
192
193//_____________________________________________________________________________
195{
196 // Destructor
197
198 delete fDirectory;
199 delete fFieldTypeCmd;
200 delete fEquationTypeCmd;
201 delete fStepperTypeCmd;
202 delete fSetMinimumStepCmd;
203 delete fSetDeltaChordCmd;
204 delete fSetDeltaOneStepCmd;
205 delete fSetDeltaIntersectionCmd;
206 delete fSetMinimumEpsilonStepCmd;
207 delete fSetMaximumEpsilonStepCmd;
208 delete fSetConstDistanceCmd;
209}
210
211//
212// public methods
213//
214
215//_____________________________________________________________________________
217 G4UIcommand* command, G4String newValues)
218{
219 // Apply command to the associated object.
220
221 if (command == fFieldTypeCmd) {
222 for (G4int i = kMagnetic; i <= kGravity; i++) {
223 G4FieldType ft = (G4FieldType)i;
224 if (newValues == G4FieldParameters::FieldTypeName(ft)) {
225 fFieldParameters->SetFieldType(ft);
226 break;
227 }
228 }
229 return;
230 }
231
232 if (command == fEquationTypeCmd) {
233 for (G4int i = kEqMagnetic; i <= kEqEMfieldWithEDM; i++) {
235 if (newValues == G4FieldParameters::EquationTypeName(et)) {
236 fFieldParameters->SetEquationType(et);
237 break;
238 }
239 }
240 return;
241 }
242
243 if (command == fStepperTypeCmd) {
244 for (G4int i = kCashKarpRKF45; i <= kRK547FEq3; i++) {
246 if (newValues == G4FieldParameters::StepperTypeName(st)) {
247 fFieldParameters->SetStepperType(st);
248 break;
249 }
250 }
251 return;
252 }
253
254 if (command == fSetMinimumStepCmd) {
255 fFieldParameters->SetMinimumStep(
256 fSetMinimumStepCmd->GetNewDoubleValue(newValues));
257 return;
258 }
259
260 if (command == fSetDeltaChordCmd) {
261 fFieldParameters->SetDeltaChord(
262 fSetDeltaChordCmd->GetNewDoubleValue(newValues));
263 return;
264 }
265
266 if (command == fSetDeltaOneStepCmd) {
267 fFieldParameters->SetDeltaOneStep(
268 fSetDeltaOneStepCmd->GetNewDoubleValue(newValues));
269 return;
270 }
271
272 if (command == fSetDeltaIntersectionCmd) {
273 fFieldParameters->SetDeltaIntersection(
274 fSetDeltaIntersectionCmd->GetNewDoubleValue(newValues));
275 return;
276 }
277
278 if (command == fSetMinimumEpsilonStepCmd) {
279 fFieldParameters->SetMinimumEpsilonStep(
280 fSetMinimumEpsilonStepCmd->GetNewDoubleValue(newValues));
281 return;
282 }
283
284 if (command == fSetMaximumEpsilonStepCmd) {
285 fFieldParameters->SetMaximumEpsilonStep(
286 fSetMaximumEpsilonStepCmd->GetNewDoubleValue(newValues));
287 return;
288 }
289
290 if (command == fSetConstDistanceCmd) {
291 fFieldParameters->SetConstDistance(
292 fSetConstDistanceCmd->GetNewDoubleValue(newValues));
293 return;
294 }
295
296 if (command == fPrintParametersCmd) {
297 fFieldParameters->PrintParameters();
298 return;
299 }
300}
@ G4State_Init
@ G4State_Idle
@ G4State_PreInit
Definition of the G4FieldParameters class.
G4EquationType
@ kEqMagnetic
@ kEqEMfieldWithEDM
G4FieldType
The available fields in Geant4.
@ kGravity
gravity field
@ kMagnetic
magnetic field
G4StepperType
@ kCashKarpRKF45
G4CashKarpRKF45.
@ kUserStepper
User defined stepper.
@ kRK547FEq3
G4RK547FEq3.
int G4int
Definition G4Types.hh:85
G4FieldParametersMessenger(G4FieldParameters *fieldParameters)
Standard constructor.
~G4FieldParametersMessenger() override
Destructor.
void SetNewValue(G4UIcommand *command, G4String newValues) override
Apply command to the associated object.
The magnetic field parameters.
G4String GetVolumeName() const
static G4String EquationTypeName(G4EquationType equation)
Return the equation type as a string.
static G4String StepperTypeName(G4StepperType stepper)
Return the stepper type as a string.
static G4String FieldTypeName(G4FieldType field)
Return the field type as a string.