Geant4 VMC Version 6.6
Loading...
Searching...
No Matches
TG4FieldParametersMessenger.cxx
Go to the documentation of this file.
1//------------------------------------------------
2// The Geant4 Virtual Monte Carlo package
3// Copyright (C) 2007 - 2014 Ivana Hrivnacova
4// All rights reserved.
5//
6// For the licensing terms see geant4_vmc/LICENSE.
7// Contact: root-vmc@cern.ch
8//-------------------------------------------------
9
14
16#include "TG4FieldParameters.h"
17
18#include <G4UIcmdWithABool.hh>
19#include <G4UIcmdWithADouble.hh>
20#include <G4UIcmdWithADoubleAndUnit.hh>
21#include <G4UIcmdWithAString.hh>
22#include <G4UIcmdWithoutParameter.hh>
23#include <G4UIdirectory.hh>
24
25//_____________________________________________________________________________
27 TG4FieldParameters* fieldParameters)
28 : G4UImessenger(),
29 fFieldParameters(fieldParameters),
30 fDirectory(0),
31 fFieldTypeCmd(0),
32 fEquationTypeCmd(0),
33 fStepperTypeCmd(0),
34 fSetStepMinimumCmd(0),
35 fSetDeltaChordCmd(0),
36 fSetDeltaOneStepCmd(0),
37 fSetDeltaIntersectionCmd(0),
38 fSetMinimumEpsilonStepCmd(0),
39 fSetMaximumEpsilonStepCmd(0),
40 fSetConstDistanceCmd(0),
41 fSetIsMonopoleCmd(0),
42 fPrintParametersCmd(0)
43{
45
46 G4String directoryName = "/mcMagField/";
47 if (fieldParameters->GetVolumeName() != "") {
48 directoryName.append(fieldParameters->GetVolumeName());
49 directoryName.append("/");
50 }
51 fDirectory = new G4UIdirectory(directoryName);
52 fDirectory->SetGuidance("Magnetic field control commands.");
53
54 G4String commandName = directoryName;
55 commandName.append("fieldType");
56 fFieldTypeCmd = new G4UIcmdWithAString(commandName, this);
57 G4String guidance = "Select type of the field";
58 fFieldTypeCmd->SetGuidance(guidance);
59 fFieldTypeCmd->SetParameterName("FieldType", false);
60 G4String candidates;
61 for (G4int i = kMagnetic; i <= kGravity; i++) {
62 FieldType ft = (FieldType)i;
63 candidates += TG4FieldParameters::FieldTypeName(ft);
64 candidates += " ";
65 }
66 fFieldTypeCmd->SetCandidates(candidates);
67 fFieldTypeCmd->AvailableForStates(
68 G4State_PreInit, G4State_Init, G4State_Idle);
69
70 commandName = directoryName;
71 commandName.append("equationType");
72 fEquationTypeCmd = new G4UIcmdWithAString(commandName, this);
73 guidance = "Select type of the equation of motion of a particle in a field";
74 fEquationTypeCmd->SetGuidance(guidance);
75 fEquationTypeCmd->SetParameterName("EquationType", false);
76 candidates = "";
77 for (G4int i = kMagUsualEqRhs; i <= kEqEMFieldWithEDM; i++) {
80 candidates += " ";
81 }
82 fEquationTypeCmd->SetCandidates(candidates);
83 fEquationTypeCmd->AvailableForStates(
84 G4State_PreInit, G4State_Init, G4State_Idle);
85
86 commandName = directoryName;
87 commandName.append("stepperType");
88 fStepperTypeCmd = new G4UIcmdWithAString(commandName, this);
89 guidance =
90 "Select type of the the integrator of particle's equation of motion in a "
91 "field";
92 fStepperTypeCmd->SetGuidance(guidance);
93 fStepperTypeCmd->SetParameterName("StepperType", false);
94 candidates = "";
95 for (G4int i = kCashKarpRKF45; i <= kRKG3Stepper; i++) {
98 candidates += " ";
99 }
100 fStepperTypeCmd->SetCandidates(candidates);
101 fStepperTypeCmd->AvailableForStates(
102 G4State_PreInit, G4State_Init, G4State_Idle);
103
104 commandName = directoryName;
105 commandName.append("setStepMinimum");
106 fSetStepMinimumCmd = new G4UIcmdWithADoubleAndUnit(commandName, this);
107 fSetStepMinimumCmd->SetGuidance("Set minimum step in G4ChordFinder");
108 fSetStepMinimumCmd->SetParameterName("StepMinimum", false);
109 fSetStepMinimumCmd->SetDefaultUnit("mm");
110 fSetStepMinimumCmd->SetUnitCategory("Length");
111 fSetStepMinimumCmd->AvailableForStates(
112 G4State_PreInit, G4State_Init, G4State_Idle);
113
114 commandName = directoryName;
115 commandName.append("setDeltaChord");
116 fSetDeltaChordCmd = new G4UIcmdWithADoubleAndUnit(commandName, this);
117 fSetDeltaChordCmd->SetGuidance("Set delta chord in G4ChordFinder");
118 fSetDeltaChordCmd->SetParameterName("DeltaChord", false);
119 fSetDeltaChordCmd->SetDefaultUnit("mm");
120 fSetDeltaChordCmd->SetUnitCategory("Length");
121 fSetDeltaChordCmd->AvailableForStates(
122 G4State_PreInit, G4State_Init, G4State_Idle);
123
124 commandName = directoryName;
125 commandName.append("setDeltaOneStep");
126 fSetDeltaOneStepCmd = new G4UIcmdWithADoubleAndUnit(commandName, this);
127 fSetDeltaOneStepCmd->SetGuidance(
128 "Set delta one step in global field manager");
129 fSetDeltaOneStepCmd->SetParameterName("DeltaOneStep", false);
130 fSetDeltaOneStepCmd->SetDefaultUnit("mm");
131 fSetDeltaOneStepCmd->SetUnitCategory("Length");
132 fSetDeltaOneStepCmd->AvailableForStates(
133 G4State_PreInit, G4State_Init, G4State_Idle);
134
135 commandName = directoryName;
136 commandName.append("setDeltaIntersection");
137 fSetDeltaIntersectionCmd = new G4UIcmdWithADoubleAndUnit(commandName, this);
138 fSetDeltaIntersectionCmd->SetGuidance(
139 "Set delta intersection in global field manager");
140 fSetDeltaIntersectionCmd->SetParameterName("DeltaIntersection", false);
141 fSetDeltaIntersectionCmd->SetDefaultUnit("mm");
142 fSetDeltaIntersectionCmd->SetUnitCategory("Length");
143 fSetDeltaIntersectionCmd->AvailableForStates(
144 G4State_PreInit, G4State_Init, G4State_Idle);
145
146 commandName = directoryName;
147 commandName.append("setMinimumEpsilonStep");
148 fSetMinimumEpsilonStepCmd = new G4UIcmdWithADouble(commandName, this);
149 fSetMinimumEpsilonStepCmd->SetGuidance(
150 "Set minimum epsilon step in global field manager");
151 fSetMinimumEpsilonStepCmd->SetParameterName("MinimumEpsilonStep", false);
152 fSetMinimumEpsilonStepCmd->AvailableForStates(
153 G4State_PreInit, G4State_Init, G4State_Idle);
154
155 commandName = directoryName;
156 commandName.append("setMaximumEpsilonStep");
157 fSetMaximumEpsilonStepCmd = new G4UIcmdWithADouble(commandName, this);
158 fSetMaximumEpsilonStepCmd->SetGuidance(
159 "Set maximum epsilon step in global field manager");
160 fSetMaximumEpsilonStepCmd->SetParameterName("MaximumEpsilonStep", false);
161 fSetMaximumEpsilonStepCmd->AvailableForStates(
162 G4State_PreInit, G4State_Init, G4State_Idle);
163
164 commandName = directoryName;
165 commandName.append("setConstDistance");
166 fSetConstDistanceCmd = new G4UIcmdWithADoubleAndUnit(commandName, this);
167 fSetConstDistanceCmd->SetGuidance(
168 "Set the distance within which the field is considered constant.");
169 fSetConstDistanceCmd->SetGuidance(
170 "Non-zero value will trigger creating a cached magnetic field.");
171 fSetConstDistanceCmd->SetParameterName("ConstDistance", false);
172 fSetConstDistanceCmd->SetDefaultUnit("mm");
173 fSetConstDistanceCmd->SetUnitCategory("Length");
174 fSetConstDistanceCmd->SetRange("ConstDistance >= 0");
175 fSetConstDistanceCmd->AvailableForStates(G4State_PreInit);
176
177 commandName = directoryName;
178 commandName.append("setIsMonopole");
179 fSetIsMonopoleCmd = new G4UIcmdWithABool(commandName, this);
180 fSetIsMonopoleCmd->SetGuidance(
181 "Activate creating a special monopole field integration.");
182 fSetIsMonopoleCmd->SetParameterName("IsMonopole", false);
183 fSetIsMonopoleCmd->AvailableForStates(G4State_PreInit);
184
185 commandName = directoryName;
186 commandName.append("printParameters");
187 fPrintParametersCmd = new G4UIcmdWithoutParameter(commandName, this);
188 fPrintParametersCmd->SetGuidance("Prints all accuracy parameters.");
189 fPrintParametersCmd->AvailableForStates(
190 G4State_PreInit, G4State_Init, G4State_Idle);
191}
192
193//_____________________________________________________________________________
211
212//
213// public methods
214//
215
216//_____________________________________________________________________________
218 G4UIcommand* command, G4String newValues)
219{
221
222 if (command == fFieldTypeCmd) {
223 for (G4int i = kMagnetic; i <= kGravity; i++) {
224 FieldType ft = (FieldType)i;
225 if (newValues == TG4FieldParameters::FieldTypeName(ft)) {
227 break;
228 }
229 }
230 }
231 else if (command == fEquationTypeCmd) {
232 for (G4int i = kMagUsualEqRhs; i <= kEqEMFieldWithEDM; i++) {
234 if (newValues == TG4FieldParameters::EquationTypeName(et)) {
236 break;
237 }
238 }
239 }
240 else if (command == fStepperTypeCmd) {
241 for (G4int i = kCashKarpRKF45; i <= kRKG3Stepper; i++) {
242 StepperType st = (StepperType)i;
243 if (newValues == TG4FieldParameters::StepperTypeName(st)) {
245 break;
246 }
247 }
248 }
249 else if (command == fSetStepMinimumCmd) {
251 fSetStepMinimumCmd->GetNewDoubleValue(newValues));
252 }
253 else if (command == fSetDeltaChordCmd) {
255 fSetDeltaChordCmd->GetNewDoubleValue(newValues));
256 }
257 else if (command == fSetDeltaOneStepCmd) {
259 fSetDeltaOneStepCmd->GetNewDoubleValue(newValues));
260 }
261 else if (command == fSetDeltaIntersectionCmd) {
263 fSetDeltaIntersectionCmd->GetNewDoubleValue(newValues));
264 }
265 else if (command == fSetMinimumEpsilonStepCmd) {
267 fSetMinimumEpsilonStepCmd->GetNewDoubleValue(newValues));
268 }
269 else if (command == fSetMaximumEpsilonStepCmd) {
271 fSetMaximumEpsilonStepCmd->GetNewDoubleValue(newValues));
272 }
273 else if (command == fSetConstDistanceCmd) {
275 fSetConstDistanceCmd->GetNewDoubleValue(newValues));
276 }
277 else if (command == fSetIsMonopoleCmd) {
279 fSetIsMonopoleCmd->GetNewBoolValue(newValues));
280 }
281 else if (command == fPrintParametersCmd) {
283 }
284}
Definition of the TG4FieldParametersMessenger class.
Definition of the TG4FieldParameters class.
@ kRKG3Stepper
G4RKG3_Stepper.
@ kCashKarpRKF45
G4CashKarpRKF45.
@ kMagUsualEqRhs
@ kEqEMFieldWithEDM
FieldType
The available fields in Geant4.
@ kGravity
gravity field
@ kMagnetic
magnetic field
G4UIcmdWithoutParameter * fPrintParametersCmd
command: printParameters
G4UIcmdWithADouble * fSetMinimumEpsilonStepCmd
command: setMinimumEpsilon
G4UIcmdWithADoubleAndUnit * fSetDeltaChordCmd
command: setDeltaChord
G4UIcmdWithADoubleAndUnit * fSetConstDistanceCmd
command: setConstDistance
G4UIcmdWithADouble * fSetMaximumEpsilonStepCmd
command: setMaximumEpsilon
G4UIcmdWithAString * fEquationTypeCmd
command: equationType
G4UIdirectory * fDirectory
command directory
virtual void SetNewValue(G4UIcommand *command, G4String newValues)
G4UIcmdWithAString * fStepperTypeCmd
command: stepperType
TG4FieldParameters * fFieldParameters
associated class
G4UIcmdWithADoubleAndUnit * fSetDeltaIntersectionCmd
command: setDeltaIntersection
G4UIcmdWithADoubleAndUnit * fSetStepMinimumCmd
command: setStepMinimum
G4UIcmdWithADoubleAndUnit * fSetDeltaOneStepCmd
command: setDeltaOneStep
G4UIcmdWithABool * fSetIsMonopoleCmd
command: setIsMonopole
TG4FieldParametersMessenger()
Not implemented.
G4UIcmdWithAString * fFieldTypeCmd
command: fieldType
The magnetic field parameters.
void SetEquationType(EquationType equation)
Set the type of equation of motion of a particle in a field.
void SetConstDistance(G4double value)
Set the distance within which the field is considered constant.
void SetStepperType(StepperType stepper)
Set the type of integrator of particle's equation of motion.
void SetIsMonopole(G4bool isMonopole)
G4String GetVolumeName() const
Return the name of associated volume, if local field.
void SetMaximumEpsilonStep(G4double value)
Set maximum epsilon step in global field manager.
static G4String EquationTypeName(EquationType equation)
void SetStepMinimum(G4double value)
Set minimum step in G4ChordFinder.
void SetDeltaChord(G4double value)
Set delta chord in G4ChordFinder.
static G4String StepperTypeName(StepperType stepper)
void SetDeltaIntersection(G4double value)
Set delta intersection in global field manager.
void SetFieldType(FieldType field)
Set type of field.
void SetDeltaOneStep(G4double value)
Set delta one step in global field manager.
void SetMinimumEpsilonStep(G4double value)
Set minimum epsilon step in global field manager.
static G4String FieldTypeName(FieldType field)