Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4FieldParametersMessenger.hh
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/// \file G4FieldParametersMessenger.h
27/// \brief Definition of the G4FieldParametersMessenger class
28///
29/// This code was initially developed in Geant4 VMC package
30/// (https://github.com/vmc-project)
31/// and adapted to Geant4.
32///
33/// \author I. Hrivnacova; IJCLab, Orsay
34
35#ifndef G4FIELDPARAMETERSMESSENGER_HH
36#define G4FIELDPARAMETERSMESSENGER_HH
37
38#include "G4UImessenger.hh"
39#include "globals.hh"
40
42
43class G4UIcommand;
44class G4UIdirectory;
50
51/// \ingroup geometry
52/// \brief Messenger class that defines commands for TG4DetConstruction.
53///
54/// Implements commands:
55/// - /field/fieldType fieldType \n
56/// fieldType = Magnetic | ElectroMagnetic | Gravity
57/// - /field/equationType eqType \n
58/// eqType = EqMagnetic | EqMagneticWithSpin | EqElectroMagnetic |
59/// EqEMfieldWithSpin | EqEMfieldWithEDM
60/// - /field/stepperType stepperType \n
61/// stepperType = CashKarpRKF45 | ClassicalRK4 | ExplicitEuler | ImplicitEuler |
62/// SimpleHeum | SimpleRunge | ConstRK4 | ExactHelixStepper
63/// | HelixExplicitEuler | HelixHeum | HelixImplicitEuler |
64/// HelixMixedStepper | HelixSimpleRunge | NystromRK4 |
65/// RKG3Stepper
66/// - /field/setMinimumStep value
67/// - /field/setDeltaChord value
68/// - /field/setDeltaOneStep value
69/// - /field/setDeltaIntersection value
70/// - /field/setMinimumEpsilonStep value
71/// - /field/setMaximumEpsilonStep value
72/// - /field/setConstDistance value
73/// - /field/printParameters
74///
75/// \author I. Hrivnacova; IJClab, Orsay
76
78{
79 public:
80 /// Standard constructor
82 /// Destructor
84
85 // methods
86 /// Apply command to the associated object.
87 void SetNewValue(G4UIcommand* command, G4String newValues) override;
88
89 private:
90 /// Not implemented
92 /// Not implemented
94 /// Not implemented
96 const G4FieldParametersMessenger& right) = delete;
97
98 // Data members
99
100 G4FieldParameters* fFieldParameters = nullptr; ///< associated class
101 G4UIdirectory* fDirectory = nullptr; ///< command directory
102
103 // Commands data members
104
105 /// Command: fieldType
106 G4UIcmdWithAString* fFieldTypeCmd = nullptr;
107
108 /// Command: equationType
109 G4UIcmdWithAString* fEquationTypeCmd = nullptr;
110
111 /// Command: stepperType
112 G4UIcmdWithAString* fStepperTypeCmd = nullptr;
113
114 /// Command: setMinimumStep
115 G4UIcmdWithADoubleAndUnit* fSetMinimumStepCmd = nullptr;
116
117 /// Command: setDeltaChord
118 G4UIcmdWithADoubleAndUnit* fSetDeltaChordCmd = nullptr;
119
120 /// Command: setDeltaOneStep
121 G4UIcmdWithADoubleAndUnit* fSetDeltaOneStepCmd = nullptr;
122
123 /// Command: setDeltaIntersection
124 G4UIcmdWithADoubleAndUnit* fSetDeltaIntersectionCmd = nullptr;
125
126 /// Command: setMinimumEpsilon
127 G4UIcmdWithADouble* fSetMinimumEpsilonStepCmd = nullptr;
128
129 /// Command: setMaximumEpsilon
130 G4UIcmdWithADouble* fSetMaximumEpsilonStepCmd = nullptr;
131
132 /// Command: setConstDistance
133 G4UIcmdWithADoubleAndUnit* fSetConstDistanceCmd = nullptr;
134
135 /// Command: printParameters
136 G4UIcmdWithoutParameter* fPrintParametersCmd = nullptr;
137};
138
139#endif // G4FIELDPARAMETERSMESSENGER_HH
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.
G4UImessenger()=default