CMS 3D CMS Logo

CMSFieldManager.cc
Go to the documentation of this file.
4 
5 #include "CLHEP/Units/GlobalSystemOfUnits.h"
6 #include "G4ChordFinder.hh"
7 #include "G4MagIntegratorStepper.hh"
8 #include "G4PropagatorInField.hh"
9 #include "G4Region.hh"
10 #include "G4RegionStore.hh"
11 #include "G4Track.hh"
12 
14  : G4FieldManager(),
15  m_currChordFinder(nullptr),
16  m_chordFinder(nullptr),
17  m_chordFinderMonopole(nullptr),
18  m_propagator(nullptr),
19  m_dChord(0.001),
20  m_dChordTracker(0.001),
21  m_dOneStep(0.001),
22  m_dOneStepTracker(0.0001),
23  m_dIntersection(0.0001),
24  m_dInterTracker(1e-6),
25  m_Rmax2(1.e+6),
26  m_Zmax(3.e+3),
27  m_stepMax(1000000.),
28  m_energyThTracker(1.e+7),
29  m_energyThreshold(0.0),
30  m_dChordSimple(0.1),
31  m_dOneStepSimple(0.1),
32  m_dIntersectionSimple(0.01),
33  m_stepMaxSimple(1000.),
34  m_cfTracker(false),
35  m_cfVacuum(false) {}
36 
39  delete m_chordFinder;
40  }
42  delete m_chordFinderMonopole;
43  }
44 }
45 
47  sim::Field *field,
48  G4ChordFinder *cf,
49  G4ChordFinder *cfmon,
50  const std::string &vol,
51  const std::string &type,
52  const std::string &stepper,
53  double delta,
54  G4PropagatorInField *pf) {
55  double minstep = p.getParameter<double>("MinStep") * CLHEP::mm;
56  double minEpsStep = p.getUntrackedParameter<double>("MinimumEpsilonStep", 0.00001) * CLHEP::mm;
57  double maxEpsStep = p.getUntrackedParameter<double>("MaximumEpsilonStep", 0.01) * CLHEP::mm;
58  int maxLC = (int)p.getUntrackedParameter<double>("MaximumLoopCounts", 1000.);
59 
60  // double
61  m_dChord = p.getParameter<double>("DeltaChord") * CLHEP::mm;
62  m_dChordTracker = p.getParameter<double>("DeltaChord") * CLHEP::mm;
63  m_dOneStep = p.getParameter<double>("DeltaOneStep") * CLHEP::mm;
64  m_dOneStepTracker = p.getParameter<double>("DeltaOneStepTracker") * CLHEP::mm;
65  m_dIntersection = p.getParameter<double>("DeltaIntersection") * CLHEP::mm;
66  m_dInterTracker = p.getParameter<double>("DeltaIntersectionTracker") * CLHEP::mm;
67  m_stepMax = p.getParameter<double>("MaxStep") * CLHEP::cm;
68 
69  m_energyThreshold = p.getParameter<double>("EnergyThSimple") * CLHEP::GeV;
70  m_energyThTracker = p.getParameter<double>("EnergyThTracker") * CLHEP::GeV;
71 
72  double rmax = p.getParameter<double>("RmaxTracker") * CLHEP::mm;
73  m_Rmax2 = rmax * rmax;
74  m_Zmax = p.getParameter<double>("ZmaxTracker") * CLHEP::mm;
75 
76  m_dChordSimple = p.getParameter<double>("DeltaChordSimple") * CLHEP::mm;
77  m_dOneStepSimple = p.getParameter<double>("DeltaOneStepSimple") * CLHEP::mm;
78  m_dIntersectionSimple = p.getParameter<double>("DeltaIntersectionSimple") * CLHEP::mm;
79  m_stepMaxSimple = p.getParameter<double>("MaxStepSimple") * CLHEP::cm;
80 
81  edm::LogVerbatim("SimG4CoreApplication")
82  << " New CMSFieldManager: LogicalVolume: <" << vol << ">\n"
83  << " Stepper: <" << stepper << ">\n"
84  << " Field type <" << type << ">\n"
85  << " Field const delta " << delta << " mm\n"
86  << " MaximumLoopCounts " << maxLC << "\n"
87  << " MinimumEpsilonStep " << minEpsStep << "\n"
88  << " MaximumEpsilonStep " << maxEpsStep << "\n"
89  << " MinStep " << minstep << " mm\n"
90  << " MaxStep " << m_stepMax / CLHEP::cm << " cm\n"
91  << " DeltaChord " << m_dChord << " mm\n"
92  << " DeltaOneStep " << m_dOneStep << " mm\n"
93  << " DeltaIntersection " << m_dIntersection << " mm\n"
94  << " DeltaInterTracker " << m_dInterTracker << " mm\n"
95  << " EnergyThresholdSimple " << m_energyThreshold / CLHEP::MeV << " MeV\n"
96  << " EnergyThresholdTracker " << m_energyThTracker / CLHEP::MeV << " MeV\n"
97  << " DeltaChordSimple " << m_dChordSimple << " mm\n"
98  << " DeltaOneStepSimple " << m_dOneStepSimple << " mm\n"
99  << " DeltaIntersectionSimple " << m_dIntersectionSimple << " mm\n"
100  << " MaxStepInVacuum " << m_stepMaxSimple / CLHEP::cm << " cm";
101 
102  // initialisation of chord finders
103  m_chordFinder = cf;
104  m_chordFinderMonopole = cfmon;
105 
106  m_chordFinderMonopole->SetDeltaChord(m_dChord);
107 
108  // initialisation of field manager
109  theField.reset(field);
110  SetDetectorField(field);
111  SetMinimumEpsilonStep(minEpsStep);
112  SetMaximumEpsilonStep(maxEpsStep);
113 
114  // propagater in field
115  m_propagator = pf;
116  pf->SetMaxLoopCount(maxLC);
117  pf->SetMinimumEpsilonStep(minEpsStep);
118  pf->SetMaximumEpsilonStep(maxEpsStep);
119 
120  // initial initialisation the default chord finder
121  setMonopoleTracking(false);
122 
123  // define regions
124  std::vector<std::string> rnames = p.getParameter<std::vector<std::string>>("VacRegions");
125  if (!rnames.empty()) {
126  std::stringstream ss;
127  std::vector<G4Region *> *rs = G4RegionStore::GetInstance();
128  for (auto &regnam : rnames) {
129  for (auto &reg : *rs) {
130  if (regnam == reg->GetName()) {
131  m_regions.push_back(reg);
132  ss << " " << regnam;
133  }
134  }
135  }
136  edm::LogVerbatim("SimG4CoreApplication") << "Simple field integration in G4Regions:\n" << ss.str() << "\n";
137  }
138 }
139 
141  // run time parameters per track
142  if (track->GetKineticEnergy() > m_energyThTracker && isInsideTracker(track)) {
143  if (!m_cfTracker) {
145  }
146 
147  } else if (track->GetKineticEnergy() <= m_energyThreshold || isInsideVacuum(track)) {
148  if (!m_cfVacuum) {
150  }
151 
152  } else if (m_cfTracker || m_cfVacuum) {
153  // restore defaults
155  }
156 }
157 
159  if (flag) {
161  if (m_cfTracker || m_cfVacuum) {
163  }
165  SetChordFinder(m_currChordFinder);
166  }
167  } else {
169  }
170  SetFieldChangesEnergy(flag);
171  m_currChordFinder->ResetStepEstimate();
172 }
173 
175  if (!m_regions.empty()) {
176  const G4Region *reg = track->GetVolume()->GetLogicalVolume()->GetRegion();
177  for (auto &areg : m_regions) {
178  if (reg == areg) {
179  return true;
180  }
181  }
182  }
183  return false;
184 }
185 
187  const G4ThreeVector &pos = track->GetPosition();
188  const double x = pos.x();
189  const double y = pos.y();
190  return (x * x + y * y < m_Rmax2 && std::abs(pos.z()) < m_Zmax);
191 }
192 
196  SetChordFinder(m_currChordFinder);
197  }
198  m_currChordFinder->SetDeltaChord(m_dChord);
199  SetDeltaOneStep(m_dOneStep);
200  SetDeltaIntersection(m_dIntersection);
201  m_propagator->SetLargestAcceptableStep(m_stepMax);
202  m_cfVacuum = m_cfTracker = false;
203 }
204 
208  SetChordFinder(m_currChordFinder);
209  }
210  m_currChordFinder->SetDeltaChord(m_dChordTracker);
211  SetDeltaOneStep(m_dOneStepTracker);
212  SetDeltaIntersection(m_dInterTracker);
213  m_propagator->SetLargestAcceptableStep(m_stepMax);
214  m_cfVacuum = false;
215  m_cfTracker = true;
216 }
217 
221  SetChordFinder(m_currChordFinder);
222  }
223  m_currChordFinder->SetDeltaChord(m_dChordSimple);
224  SetDeltaOneStep(m_dOneStepSimple);
225  SetDeltaIntersection(m_dIntersectionSimple);
226  m_propagator->SetLargestAcceptableStep(m_stepMaxSimple);
227  m_cfVacuum = true;
228  m_cfTracker = false;
229 }
CMSFieldManager::ConfigureForTrack
void ConfigureForTrack(const G4Track *) override
Definition: CMSFieldManager.cc:140
CMSFieldManager::m_dInterTracker
double m_dInterTracker
Definition: CMSFieldManager.h:69
DDAxes::y
CMSFieldManager::m_dOneStepSimple
double m_dOneStepSimple
Definition: CMSFieldManager.h:76
CMSFieldManager::m_dOneStep
double m_dOneStep
Definition: CMSFieldManager.h:66
HLT_FULL_cff.track
track
Definition: HLT_FULL_cff.py:11713
MessageLogger.h
funct::false
false
Definition: Factorize.h:29
CMSFieldManager::m_chordFinder
G4ChordFinder * m_chordFinder
Definition: CMSFieldManager.h:57
CMSFieldManager::m_dIntersectionSimple
double m_dIntersectionSimple
Definition: CMSFieldManager.h:77
CMSFieldManager::m_stepMax
double m_stepMax
Definition: CMSFieldManager.h:72
CMSFieldManager::m_stepMaxSimple
double m_stepMaxSimple
Definition: CMSFieldManager.h:78
CMSFieldManager.h
CMSFieldManager::m_Rmax2
double m_Rmax2
Definition: CMSFieldManager.h:70
CMSFieldManager::InitialiseForVolume
void InitialiseForVolume(const edm::ParameterSet &, sim::Field *, G4ChordFinder *cfDefault, G4ChordFinder *cfMonopole, const std::string &vol, const std::string &fieldType, const std::string &stepperName, double delta, G4PropagatorInField *)
Definition: CMSFieldManager.cc:46
AlCaHLTBitMon_ParallelJobs.p
p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
pos
Definition: PixelAliasList.h:18
MeV
const double MeV
DDAxes::x
CMSFieldManager::isInsideTracker
bool isInsideTracker(const G4Track *)
Definition: CMSFieldManager.cc:186
CMSFieldManager::theField
std::unique_ptr< sim::Field > theField
Definition: CMSFieldManager.h:54
CMSFieldManager::isInsideVacuum
bool isInsideVacuum(const G4Track *)
Definition: CMSFieldManager.cc:174
CMSFieldManager::m_Zmax
double m_Zmax
Definition: CMSFieldManager.h:71
contentValuesCheck.ss
ss
Definition: contentValuesCheck.py:33
CMSFieldManager::setMonopoleTracking
void setMonopoleTracking(G4bool)
Definition: CMSFieldManager.cc:158
CMSFieldManager::m_currChordFinder
G4ChordFinder * m_currChordFinder
Definition: CMSFieldManager.h:56
CMSFieldManager::m_chordFinderMonopole
G4ChordFinder * m_chordFinderMonopole
Definition: CMSFieldManager.h:58
sim::Field
Definition: Field.h:9
CMSFieldManager::m_dChordTracker
double m_dChordTracker
Definition: CMSFieldManager.h:65
CMSFieldManager::m_dChord
double m_dChord
Definition: CMSFieldManager.h:64
CMSFieldManager::m_cfTracker
bool m_cfTracker
Definition: CMSFieldManager.h:80
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
CMSFieldManager::~CMSFieldManager
~CMSFieldManager() override
Definition: CMSFieldManager.cc:37
CMSFieldManager::m_propagator
G4PropagatorInField * m_propagator
Definition: CMSFieldManager.h:60
edm::ParameterSet
Definition: ParameterSet.h:47
CMSFieldManager::setChordFinderForTracker
void setChordFinderForTracker()
Definition: CMSFieldManager.cc:205
type
type
Definition: SiPixelVCal_PayloadInspector.cc:37
GeV
const double GeV
Definition: MathUtil.h:16
dumpMFGeometry_cfg.delta
delta
Definition: dumpMFGeometry_cfg.py:25
Field.h
createfilelist.int
int
Definition: createfilelist.py:10
CMSFieldManager::setChordFinderForVacuum
void setChordFinderForVacuum()
Definition: CMSFieldManager.cc:218
CMSFieldManager::setDefaultChordFinder
void setDefaultChordFinder()
Definition: CMSFieldManager.cc:193
packedPFCandidateRefMixer_cfi.pf
pf
Definition: packedPFCandidateRefMixer_cfi.py:4
edm::LogVerbatim
Log< level::Info, true > LogVerbatim
Definition: MessageLogger.h:128
CMSFieldManager::m_energyThTracker
double m_energyThTracker
Definition: CMSFieldManager.h:73
CMSFieldManager::m_dChordSimple
double m_dChordSimple
Definition: CMSFieldManager.h:75
CMSFieldManager::m_dOneStepTracker
double m_dOneStepTracker
Definition: CMSFieldManager.h:67
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
CMSFieldManager::CMSFieldManager
CMSFieldManager()
Definition: CMSFieldManager.cc:13
CMSFieldManager::m_regions
std::vector< const G4Region * > m_regions
Definition: CMSFieldManager.h:62
CMSFieldManager::m_energyThreshold
double m_energyThreshold
Definition: CMSFieldManager.h:74
CMSFieldManager::m_dIntersection
double m_dIntersection
Definition: CMSFieldManager.h:68
RemoveAddSevLevel.flag
flag
Definition: RemoveAddSevLevel.py:116
MillePedeFileConverter_cfg.e
e
Definition: MillePedeFileConverter_cfg.py:37
CMSFieldManager::m_cfVacuum
bool m_cfVacuum
Definition: CMSFieldManager.h:81