ATLAS Offline Software
Loading...
Searching...
No Matches
G4FieldManagerToolBase Class Reference

Tool base class for implementation/retrieval of a G4FieldManager. More...

#include "G4AtlasTools/G4FieldManagerToolBase.h"

Inheritance diagram for G4FieldManagerToolBase:
Collaboration diagram for G4FieldManagerToolBase:

Public Member Functions

 G4FieldManagerToolBase (const std::string &type, const std::string &name, const IInterface *parent)
 Standard constructor.
virtual ~G4FieldManagerToolBase ()
 Destructor.
virtual StatusCode initialize () override
 Initialize method.

Protected Member Functions

G4MagIntegratorStepper * getStepper (const std::string &stepperType, G4MagneticField *field) const
 Common method to construct a stepper of requested type.
StatusCode setFieldParameters (G4FieldManager *fieldMgr) const
 Common method to apply configuredfield parameters.

Protected Attributes

ServiceHandle< IG4FieldSvcm_fieldSvc {this, "FieldSvc", "G4FieldSvc", "Service providing a G4MagneticField"}
 Handle to the G4 field service.
Gaudi::Property< bool > m_fieldOn {this, "FieldOn", true, "Toggles field on/off"}
Gaudi::Property< std::string > m_integratorStepper {this, "IntegratorStepper", "AtlasRK4", "Integrator stepper name"}
 The type of stepper to use.
ToolHandle< IEquationOfMotionToolm_equationOfMotion {this, "EquationOfMotion", "", ""}
 The type of equation of motion to use.
Field parameters
Gaudi::Property< double > m_minEps {this, "MinimumEpsilonStep", -1.0, "Minimum epsilon (see G4 documentation)"}
Gaudi::Property< double > m_maxEps {this, "MaximumEpsilonStep", -1.0, "Maximum epsilon (see G4 documentation)"}
Gaudi::Property< double > m_deltaChord {this, "DeltaChord", -1.0, "Missing distance for the chord finder"}
Gaudi::Property< double > m_deltaOneStep {this, "DeltaOneStep", -1.0, "Delta(one-step)"}
Gaudi::Property< double > m_deltaIntersection {this, "DeltaIntersection", -1.0, "Accuracy for boundary intersection"}
Gaudi::Property< double > m_maxStep {this, "MaximumStep", -1.0, "Maximum step length in field (see G4 documentation)"}
Gaudi::Property< double > m_minStep {this, "MinimumStep",1e-2, "Minimum step length in field (see G4 documentation)"}

Detailed Description

Tool base class for implementation/retrieval of a G4FieldManager.

Author
Andrea Dell'Acqua
Date
2015-11-17

Definition at line 33 of file G4FieldManagerToolBase.h.

Constructor & Destructor Documentation

◆ G4FieldManagerToolBase()

G4FieldManagerToolBase::G4FieldManagerToolBase ( const std::string & type,
const std::string & name,
const IInterface * parent )

Standard constructor.

Definition at line 52 of file G4FieldManagerToolBase.cxx.

55 : base_class(type, name, parent)
56{
57}

◆ ~G4FieldManagerToolBase()

virtual G4FieldManagerToolBase::~G4FieldManagerToolBase ( )
inlinevirtual

Destructor.

Definition at line 42 of file G4FieldManagerToolBase.h.

42{}

Member Function Documentation

◆ getStepper()

G4MagIntegratorStepper * G4FieldManagerToolBase::getStepper ( const std::string & stepperType,
G4MagneticField * field ) const
protected

Common method to construct a stepper of requested type.

Definition at line 190 of file G4FieldManagerToolBase.cxx.

191{
192 ATH_MSG_DEBUG("getStepper");
193 G4Mag_EqRhs* eqRhs(nullptr);
194 if (!m_equationOfMotion.empty())
195 {
196 eqRhs = m_equationOfMotion->makeEquationOfMotion(field);
197 ATH_MSG_INFO("Configuring alternative equation of motion using " <<
198 m_equationOfMotion.name() );
199 }
200 else
201 {
202 ATH_MSG_VERBOSE("Using G4Mag_UsualEqRhs as the equation of motion.");
203 eqRhs = new G4Mag_UsualEqRhs(field);
204 }
205 // @TODO Add some way of confirming that the choices of Equation of
206 // motion and stepper are compatible.
207 // @TODO consider moving the stepper choice into a tool as well?
208 // Only worthwhile if we envisage adding more stepper choices in the
209 // future.
210 if (name=="HelixImplicitEuler") return new G4HelixImplicitEuler(eqRhs);
211 else if (name=="HelixSimpleRunge") return new G4HelixSimpleRunge(eqRhs);
212 else if (name=="HelixExplicitEuler") return new G4HelixExplicitEuler(eqRhs);
213 else if (name=="NystromRK4") return new G4NystromRK4(eqRhs);
214 else if (name=="ClassicalRK4") return new G4ClassicalRK4(eqRhs);
215 else if (name=="AtlasRK4") return new G4AtlasRK4(eqRhs);
216 else if (name=="BogackiShampine23") return new G4BogackiShampine23(eqRhs);
217 else if (name=="BogackiShampine45") return new G4BogackiShampine45(eqRhs);
218 else if (name=="CashKarpRKF45") return new G4CashKarpRKF45(eqRhs);
219 else if (name=="DoLoMcPriRK34") return new G4DoLoMcPriRK34(eqRhs);
220 else if (name=="DormandPrince745") return new G4DormandPrince745(eqRhs);
221 else if (name=="DormandPrinceRK56") return new G4DormandPrinceRK56(eqRhs);
222 else if (name=="DormandPrinceRK78") return new G4DormandPrinceRK78(eqRhs);
223 else if (name=="RK547FEq1") return new G4RK547FEq1(eqRhs);
224 else if (name=="RK547FEq2") return new G4RK547FEq2(eqRhs);
225 else if (name=="RK547FEq3") return new G4RK547FEq3(eqRhs);
226 else if (name=="RKG3_Stepper") return new G4RKG3_Stepper(eqRhs);
227 else if (name=="TsitourasRK45") return new G4TsitourasRK45(eqRhs);
228 else {
229 ATH_MSG_ERROR("Stepper " << name << " not available! returning NystromRK4!");
230 return new G4NystromRK4(eqRhs);
231 }
232}
#define ATH_MSG_ERROR(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_DEBUG(x)
ToolHandle< IEquationOfMotionTool > m_equationOfMotion
The type of equation of motion to use.

◆ initialize()

StatusCode G4FieldManagerToolBase::initialize ( )
overridevirtual

Initialize method.

Definition at line 62 of file G4FieldManagerToolBase.cxx.

63{
64 ATH_MSG_DEBUG("G4FieldManagerToolBase::initialize()");
65
66 // Skip field svc retrieval if field disabled.
67 if(m_fieldOn) {
68 ATH_CHECK( m_fieldSvc.retrieve() );
69 }
70 else {
71 ATH_MSG_DEBUG("Field disabled. Not retrieving G4 field svc.");
72 }
73 if (!m_equationOfMotion.empty())
74 {
75 ATH_CHECK(m_equationOfMotion.retrieve());
76 }
77 return StatusCode::SUCCESS;
78}
#define ATH_CHECK
Evaluate an expression and check for errors.
ServiceHandle< IG4FieldSvc > m_fieldSvc
Handle to the G4 field service.
Gaudi::Property< bool > m_fieldOn

◆ setFieldParameters()

StatusCode G4FieldManagerToolBase::setFieldParameters ( G4FieldManager * fieldMgr) const
protected

Common method to apply configuredfield parameters.

Definition at line 238 of file G4FieldManagerToolBase.cxx.

240{
241 ATH_MSG_DEBUG("setFieldParameters");
242 if (fieldMgr) {
243 if (m_deltaChord>0) fieldMgr->GetChordFinder()->SetDeltaChord(m_deltaChord);
244 if (m_deltaOneStep>0) fieldMgr->SetDeltaOneStep(m_deltaOneStep);
245 if (m_deltaIntersection>0) fieldMgr->SetDeltaIntersection(m_deltaIntersection);
246
247 auto minEps_actual = m_minEps > 0 ? m_minEps.value() : fieldMgr->GetMinimumEpsilonStep();
248 auto maxEps_actual = m_maxEps > 0 ? m_maxEps.value() : fieldMgr->GetMaximumEpsilonStep();
249 G4FieldManagerHelper::SetMinAndMaxEpsilonStep(fieldMgr, minEps_actual, maxEps_actual);
250 }
251 else {
252 ATH_MSG_ERROR("setFieldParameters received NULL field mgr!");
253 return StatusCode::FAILURE;
254 }
255 return StatusCode::SUCCESS;
256}
Gaudi::Property< double > m_deltaChord
Gaudi::Property< double > m_maxEps
Gaudi::Property< double > m_deltaIntersection
Gaudi::Property< double > m_deltaOneStep
Gaudi::Property< double > m_minEps
G4bool SetMinAndMaxEpsilonStep(G4FieldManager *fieldMgr, double eps_min, double eps_max)
Set epsilon step range for a G4FieldManager instance accounting for Geant4 sanity checks.

Member Data Documentation

◆ m_deltaChord

Gaudi::Property<double> G4FieldManagerToolBase::m_deltaChord {this, "DeltaChord", -1.0, "Missing distance for the chord finder"}
protected

Definition at line 81 of file G4FieldManagerToolBase.h.

81{this, "DeltaChord", -1.0, "Missing distance for the chord finder"};

◆ m_deltaIntersection

Gaudi::Property<double> G4FieldManagerToolBase::m_deltaIntersection {this, "DeltaIntersection", -1.0, "Accuracy for boundary intersection"}
protected

Definition at line 83 of file G4FieldManagerToolBase.h.

83{this, "DeltaIntersection", -1.0, "Accuracy for boundary intersection"};

◆ m_deltaOneStep

Gaudi::Property<double> G4FieldManagerToolBase::m_deltaOneStep {this, "DeltaOneStep", -1.0, "Delta(one-step)"}
protected

Definition at line 82 of file G4FieldManagerToolBase.h.

82{this, "DeltaOneStep", -1.0, "Delta(one-step)"};

◆ m_equationOfMotion

ToolHandle<IEquationOfMotionTool> G4FieldManagerToolBase::m_equationOfMotion {this, "EquationOfMotion", "", ""}
protected

The type of equation of motion to use.

Definition at line 74 of file G4FieldManagerToolBase.h.

74{this, "EquationOfMotion", "", ""};

◆ m_fieldOn

Gaudi::Property<bool> G4FieldManagerToolBase::m_fieldOn {this, "FieldOn", true, "Toggles field on/off"}
protected
Todo
TODO why is this duplicated in the g4 field svc?

Definition at line 68 of file G4FieldManagerToolBase.h.

68{this, "FieldOn", true, "Toggles field on/off"};

◆ m_fieldSvc

ServiceHandle<IG4FieldSvc> G4FieldManagerToolBase::m_fieldSvc {this, "FieldSvc", "G4FieldSvc", "Service providing a G4MagneticField"}
protected

Handle to the G4 field service.

Definition at line 65 of file G4FieldManagerToolBase.h.

65{this, "FieldSvc", "G4FieldSvc", "Service providing a G4MagneticField"};

◆ m_integratorStepper

Gaudi::Property<std::string> G4FieldManagerToolBase::m_integratorStepper {this, "IntegratorStepper", "AtlasRK4", "Integrator stepper name"}
protected

The type of stepper to use.

Definition at line 71 of file G4FieldManagerToolBase.h.

71{this, "IntegratorStepper", "AtlasRK4", "Integrator stepper name"};

◆ m_maxEps

Gaudi::Property<double> G4FieldManagerToolBase::m_maxEps {this, "MaximumEpsilonStep", -1.0, "Maximum epsilon (see G4 documentation)"}
protected

Definition at line 80 of file G4FieldManagerToolBase.h.

80{this, "MaximumEpsilonStep", -1.0, "Maximum epsilon (see G4 documentation)"};

◆ m_maxStep

Gaudi::Property<double> G4FieldManagerToolBase::m_maxStep {this, "MaximumStep", -1.0, "Maximum step length in field (see G4 documentation)"}
protected

Definition at line 84 of file G4FieldManagerToolBase.h.

84{this, "MaximumStep", -1.0, "Maximum step length in field (see G4 documentation)"};

◆ m_minEps

Gaudi::Property<double> G4FieldManagerToolBase::m_minEps {this, "MinimumEpsilonStep", -1.0, "Minimum epsilon (see G4 documentation)"}
protected

Definition at line 79 of file G4FieldManagerToolBase.h.

79{this, "MinimumEpsilonStep", -1.0, "Minimum epsilon (see G4 documentation)"};

◆ m_minStep

Gaudi::Property<double> G4FieldManagerToolBase::m_minStep {this, "MinimumStep",1e-2, "Minimum step length in field (see G4 documentation)"}
protected

Definition at line 85 of file G4FieldManagerToolBase.h.

85{this, "MinimumStep",1e-2, "Minimum step length in field (see G4 documentation)"};

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