ATLAS Offline Software
Loading...
Searching...
No Matches
SwitchingFieldManager.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2021 CERN for the benefit of the ATLAS collaboration
3*/
4
5#include "G4ThreeVector.hh"
6#include "G4Field.hh"
7
10#include "G4ChordFinder.hh"
11
12#include "G4Track.hh"
13#include "G4MuonPlus.hh"
14#include "G4MuonMinus.hh"
15#include "G4ChargedGeantino.hh"
16
17
18
20 : G4FieldManager( detectorField ),
21 m_savedField( detectorField )
22{}
23
26
27// Inside Calo area switch off the B-field.
28//
30{
31 assert( GetDetectorField() == m_savedField
32 || GetDetectorField() == nullptr );
33 // Must ensure that SetDetectorField() has NOT been called to establish
34 // a different field by an inside entity
35 // ... it would invalidate our assumptions !
36
37 const bool isMuonTrack = (track->GetDefinition()==G4MuonPlus::Definition() ||
38 track->GetDefinition()==G4MuonMinus::Definition() ||
39 track->GetDefinition()==G4ChargedGeantino::ChargedGeantinoDefinition());
40
41 const G4ThreeVector position = track->GetPosition();
42 const G4double r2XY = position.x()*position.x() + position.y()*position.y();
43 const G4double Z = position.z();
44 const G4double r2Min = (m_radiusXYmin + m_offset)*(m_radiusXYmin + m_offset);
45 const G4double r2Max = (m_radiusXYmax - m_offset)*(m_radiusXYmax - m_offset);
46
47 const bool inSide = (r2XY > r2Min) && (r2XY < r2Max) && (fabs(Z) < m_Zmax - m_offset);
48
49 if ( inSide && DoesFieldExist() && !isMuonTrack )
50 ChangeDetectorField( nullptr );
51 else if ( !inSide && !DoesFieldExist() )
52 ChangeDetectorField( m_savedField );
53}
54
55G4FieldManager* SwitchingFieldManager::Clone() const
56{
57 G4Field* cloneField = m_savedField->Clone();
58 auto clone= new SwitchingFieldManager( cloneField );
59
60 G4FieldManagerHelper::SetMinAndMaxEpsilonStep( clone, this->GetMinimumEpsilonStep(), this->GetMaximumEpsilonStep() );
61 clone->SetDeltaOneStep( this->GetDeltaOneStep() );
62 clone->SetDeltaIntersection( this->GetDeltaIntersection() );
63
64 clone->GetChordFinder()->SetDeltaChord( this->GetChordFinder()->GetDeltaChord() );
65 // This lives in our chord finder ...
66
67 return clone;
68}
Helper functions for G4FieldManager to adapt to changes in Geant4 interfaces since 10....
SwitchingFieldManager(G4Field *field)
G4FieldManager * Clone() const override
void ConfigureForTrack(const G4Track *) override
G4bool SetMinAndMaxEpsilonStep(G4FieldManager *fieldMgr, double eps_min, double eps_max)
Set epsilon step range for a G4FieldManager instance accounting for Geant4 sanity checks.