|
ATLAS Offline Software
|
Tool for setting up a volume-local magnetic field manager.
More...
#include "G4AtlasTools/DetectorFieldManagerTool.h"
|
Gaudi::Property< std::vector< std::string > > | m_logVolumeList {this, "LogicalVolumes", {}, "List of logical volumes to which the field will be applied"} |
| List of volumes to assign this field configuration to. More...
|
|
Gaudi::Property< std::vector< std::string > > | m_physVolumeList {this, "PhysicalVolumes", {}, "List of physical volumes to which the field will be applied"} |
|
Gaudi::Property< bool > | m_muonOnlyField {this, "MuonOnlyField", false, "Only muons experience the magnetic field"} |
| Option for muons feeling the B-field only. More...
|
|
thread_utils::ThreadLocalOwner< G4FieldManager > | m_fieldMgrHolder |
| My field manager. More...
|
|
ServiceHandle< IG4FieldSvc > | m_fieldSvc {this, "FieldSvc", "G4FieldSvc", "Service providing a G4MagneticField"} |
| Handle to the G4 field service. More...
|
|
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. More...
|
|
ToolHandle< IEquationOfMotionTool > | m_equationOfMotion {this, "EquationOfMotion", "", ""} |
| The type of equation of motion to use. More...
|
|
Tool for setting up a volume-local magnetic field manager.
- Author
- Andrea Dell'Acqua
- Date
- 2016-03-25
Definition at line 22 of file DetectorFieldManagerTool.h.
◆ DetectorFieldManagerTool()
DetectorFieldManagerTool::DetectorFieldManagerTool |
( |
const std::string & |
type, |
|
|
const std::string & |
name, |
|
|
const IInterface * |
parent |
|
) |
| |
◆ ~DetectorFieldManagerTool()
DetectorFieldManagerTool::~DetectorFieldManagerTool |
( |
| ) |
|
|
inline |
◆ getStepper()
G4MagIntegratorStepper * G4FieldManagerToolBase::getStepper |
( |
const std::string & |
stepperType, |
|
|
G4MagneticField * |
field |
|
) |
| const |
|
protectedinherited |
Common method to construct a stepper of requested type.
Definition at line 190 of file G4FieldManagerToolBase.cxx.
193 G4Mag_EqRhs* eqRhs(
nullptr);
197 ATH_MSG_INFO(
"Configuring alternative equation of motion using " <<
203 eqRhs =
new G4Mag_UsualEqRhs(
field);
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);
230 return new G4NystromRK4(eqRhs);
◆ initialize()
StatusCode G4FieldManagerToolBase::initialize |
( |
| ) |
|
|
overridevirtualinherited |
◆ initializeField()
StatusCode DetectorFieldManagerTool::initializeField |
( |
| ) |
|
|
finaloverride |
Initialize a field manager.
Definition at line 35 of file DetectorFieldManagerTool.cxx.
44 ATH_MSG_ERROR(
"DetectorFieldManagerTool::initializeField() - " <<
45 "Field manager already exists!");
46 return StatusCode::FAILURE;
53 G4FieldManager * fieldMgr =
nullptr;
58 fieldMgr =
new G4FieldManager();
65 fieldMgr->SetDetectorField(
field);
66 fieldMgr->CreateChordFinder(
field);
70 fieldMgr->CreateChordFinder(
field);
72 #if G4VERSION_NUMBER < 1040
75 G4MagInt_Driver* magDriver = fieldMgr->GetChordFinder()->GetIntegrationDriver();
76 magDriver->RenewStepperAndAdjust(stepper);
80 G4ChordFinder* chordFinder = fieldMgr->GetChordFinder();
81 chordFinder->SetIntegrationDriver(
driver);
86 auto logVolStore = G4LogicalVolumeStore::GetInstance();
88 G4LogicalVolume* logicalVolume = logVolStore->GetVolume(volume);
89 if (logicalVolume !=
nullptr) logicalVolume->SetFieldManager(fieldMgr,
true);
91 ATH_MSG_WARNING(
"No volume called " << volume <<
" was found in the G4LogicalVolumeStore! Skipping this volume.");
95 auto physVolStore = G4PhysicalVolumeStore::GetInstance();
97 G4VPhysicalVolume* physicalVolume = physVolStore->GetVolume(volume);
98 if (physicalVolume !=
nullptr) physicalVolume->GetLogicalVolume()->SetFieldManager(fieldMgr,
true);
100 ATH_MSG_WARNING(
"No volume called " << volume <<
" was found in the G4PhysicalVolumeStore! Skipping this volume.");
104 ATH_MSG_WARNING(
"No volumes are provided. Field manager is NOT assigned.");
108 return StatusCode::SUCCESS;
◆ setFieldParameters()
StatusCode G4FieldManagerToolBase::setFieldParameters |
( |
G4FieldManager * |
fieldMgr | ) |
const |
|
protectedinherited |
Common method to apply configuredfield parameters.
Definition at line 238 of file G4FieldManagerToolBase.cxx.
247 auto minEps_actual =
m_minEps > 0 ?
m_minEps.value() : fieldMgr->GetMinimumEpsilonStep();
248 auto maxEps_actual =
m_maxEps > 0 ?
m_maxEps.value() : fieldMgr->GetMaximumEpsilonStep();
252 ATH_MSG_ERROR(
"setFieldParameters received NULL field mgr!");
253 return StatusCode::FAILURE;
255 return StatusCode::SUCCESS;
◆ m_deltaChord
Gaudi::Property<double> G4FieldManagerToolBase::m_deltaChord {this, "DeltaChord", -1.0, "Missing distance for the chord finder"} |
|
protectedinherited |
◆ m_deltaIntersection
Gaudi::Property<double> G4FieldManagerToolBase::m_deltaIntersection {this, "DeltaIntersection", -1.0, "Accuracy for boundary intersection"} |
|
protectedinherited |
◆ m_deltaOneStep
Gaudi::Property<double> G4FieldManagerToolBase::m_deltaOneStep {this, "DeltaOneStep", -1.0, "Delta(one-step)"} |
|
protectedinherited |
◆ m_equationOfMotion
ToolHandle<IEquationOfMotionTool> G4FieldManagerToolBase::m_equationOfMotion {this, "EquationOfMotion", "", ""} |
|
protectedinherited |
◆ m_fieldMgrHolder
◆ m_fieldOn
Gaudi::Property<bool> G4FieldManagerToolBase::m_fieldOn {this, "FieldOn", true, "Toggles field on/off"} |
|
protectedinherited |
◆ m_fieldSvc
ServiceHandle<IG4FieldSvc> G4FieldManagerToolBase::m_fieldSvc {this, "FieldSvc", "G4FieldSvc", "Service providing a G4MagneticField"} |
|
protectedinherited |
◆ m_integratorStepper
Gaudi::Property<std::string> G4FieldManagerToolBase::m_integratorStepper {this, "IntegratorStepper", "AtlasRK4", "Integrator stepper name"} |
|
protectedinherited |
◆ m_logVolumeList
Gaudi::Property<std::vector<std::string> > DetectorFieldManagerTool::m_logVolumeList {this, "LogicalVolumes", {}, "List of logical volumes to which the field will be applied"} |
|
protected |
◆ m_maxEps
Gaudi::Property<double> G4FieldManagerToolBase::m_maxEps {this, "MaximumEpsilonStep", -1.0, "Maximum epsilon (see G4 documentation)"} |
|
protectedinherited |
◆ m_maxStep
Gaudi::Property<double> G4FieldManagerToolBase::m_maxStep {this, "MaximumStep", -1.0, "Maximum step length in field (see G4 documentation)"} |
|
protectedinherited |
◆ m_minEps
Gaudi::Property<double> G4FieldManagerToolBase::m_minEps {this, "MinimumEpsilonStep", -1.0, "Minimum epsilon (see G4 documentation)"} |
|
protectedinherited |
◆ m_minStep
Gaudi::Property<double> G4FieldManagerToolBase::m_minStep {this, "MinimumStep",1e-2, "Minimum step length in field (see G4 documentation)"} |
|
protectedinherited |
◆ m_muonOnlyField
Gaudi::Property<bool> DetectorFieldManagerTool::m_muonOnlyField {this, "MuonOnlyField", false, "Only muons experience the magnetic field"} |
|
protected |
◆ m_physVolumeList
Gaudi::Property<std::vector<std::string> > DetectorFieldManagerTool::m_physVolumeList {this, "PhysicalVolumes", {}, "List of physical volumes to which the field will be applied"} |
|
protected |
The documentation for this class was generated from the following files:
G4FieldManager that sets tight stepping for muons; disables magnetic field for other particles.
#define ATH_MSG_VERBOSE(x)
void set(T *obj)
Assign the object of the current thread.
T * get()
Get the object of the current thread.
G4bool SetMinAndMaxEpsilonStep(G4FieldManager *fieldMgr, double eps_min, double eps_max)
Set epsilon step range for a G4FieldManager instance accounting for Geant4 sanity checks.
#define ATH_MSG_WARNING(x)