ATLAS Offline Software
TightMuonSteppingFieldManager.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2017 CERN for the benefit of the ATLAS collaboration
3 */
4 
6 #include "G4Track.hh"
7 #include "G4ChordFinder.hh"
8 #include "G4MuonPlus.hh"
9 #include "G4MuonMinus.hh"
10 #include "G4ChargedGeantino.hh"
11 
12 
14 TightMuonSteppingFieldManager::TightMuonSteppingFieldManager(G4Field *detectorField, G4ChordFinder *pChordFinder, G4bool b)
15  : G4FieldManager(detectorField,pChordFinder,b)
16  , m_globalDeltaChord(-1.)
17  , m_globalDeltaOneStep(-1.)
18  , m_globalDeltaIntersection(-1.)
19  , m_globalMinEps(-1.)
20  , m_globalMaxEps(-1.)
21 {}
22 
24 TightMuonSteppingFieldManager::TightMuonSteppingFieldManager(G4MagneticField *detectorMagneticField)
25  : G4FieldManager(detectorMagneticField)
26  , m_globalDeltaChord(-1.)
27  , m_globalDeltaOneStep(-1.)
28  , m_globalDeltaIntersection(-1.)
29  , m_globalMinEps(-1.)
30  , m_globalMaxEps(-1.)
31 {}
32 
35 {
36  // If they have not been set yet, get the settings for the global field manager
37  if (m_globalDeltaChord<0){
38  m_globalDeltaChord = GetChordFinder()->GetDeltaChord();
39  m_globalDeltaOneStep = GetDeltaOneStep();
40  m_globalDeltaIntersection = GetDeltaIntersection();
41  m_globalMinEps = GetMinimumEpsilonStep();
42  m_globalMaxEps = GetMaximumEpsilonStep();
43  }
44 
45  // If this is a muon, set tight parameters; otherwise go back to the defaults
46  if (track->GetDefinition()==G4MuonPlus::Definition() ||
47  track->GetDefinition()==G4MuonMinus::Definition() ||
48  track->GetDefinition()==G4ChargedGeantino::ChargedGeantinoDefinition()) {
49  GetChordFinder()->SetDeltaChord(0.00000002);
50  SetDeltaOneStep(0.000001);
51  SetDeltaIntersection(0.00000002);
52  SetMinimumEpsilonStep(0.0000009);
53  SetMaximumEpsilonStep(0.000001);
54  } else {
55  GetChordFinder()->SetDeltaChord(m_globalDeltaChord);
56  SetDeltaOneStep(m_globalDeltaOneStep);
57  SetDeltaIntersection(m_globalDeltaIntersection);
58  SetMinimumEpsilonStep(m_globalMinEps);
59  SetMaximumEpsilonStep(m_globalMaxEps);
60  }
61 
62 }
TightMuonSteppingFieldManager.h
TightMuonSteppingFieldManager::m_globalDeltaOneStep
double m_globalDeltaOneStep
Definition: TightMuonSteppingFieldManager.h:41
TightMuonSteppingFieldManager::m_globalDeltaIntersection
double m_globalDeltaIntersection
Definition: TightMuonSteppingFieldManager.h:42
TightMuonSteppingFieldManager::m_globalMinEps
double m_globalMinEps
Definition: TightMuonSteppingFieldManager.h:43
TightMuonSteppingFieldManager::m_globalMaxEps
double m_globalMaxEps
Definition: TightMuonSteppingFieldManager.h:44
TightMuonSteppingFieldManager::m_globalDeltaChord
double m_globalDeltaChord
Parameters of the stepper.
Definition: TightMuonSteppingFieldManager.h:40
TightMuonSteppingFieldManager::TightMuonSteppingFieldManager
TightMuonSteppingFieldManager(G4Field *detectorField=0, G4ChordFinder *pChordFinder=0, G4bool b=true)
Constructor.
Definition: TightMuonSteppingFieldManager.cxx:14
plotBeamSpotMon.b
b
Definition: plotBeamSpotMon.py:77
TightMuonSteppingFieldManager::ConfigureForTrack
virtual void ConfigureForTrack(const G4Track *) override final
The one interesting method.
Definition: TightMuonSteppingFieldManager.cxx:34
xAOD::track
@ track
Definition: TrackingPrimitives.h:512