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_dOneStep(0.001),
21  m_dIntersection(0.0001),
22  m_stepMax(1000000.),
23  m_energyThreshold(0.0),
24  m_dChordSimple(0.1),
25  m_dOneStepSimple(0.1),
26  m_dIntersectionSimple(0.01),
27  m_stepMaxSimple(1000.),
28  m_cfVacuum(false) {}
29 
32  delete m_chordFinder;
33  }
35  delete m_chordFinderMonopole;
36  }
37 }
38 
40  sim::Field *field,
41  G4ChordFinder *cf,
42  G4ChordFinder *cfmon,
43  const std::string &vol,
44  const std::string &type,
45  const std::string &stepper,
46  double delta,
47  G4PropagatorInField *pf) {
48  double minstep = p.getParameter<double>("MinStep") * CLHEP::mm;
49  double minEpsStep = p.getUntrackedParameter<double>("MinimumEpsilonStep", 0.00001) * CLHEP::mm;
50  double maxEpsStep = p.getUntrackedParameter<double>("MaximumEpsilonStep", 0.01) * CLHEP::mm;
51  int maxLC = (int)p.getUntrackedParameter<double>("MaximumLoopCounts", 1000.);
52 
53  // double
54  m_dChord = p.getParameter<double>("DeltaChord") * CLHEP::mm;
55  m_dOneStep = p.getParameter<double>("DeltaOneStep") * CLHEP::mm;
56  m_dIntersection = p.getParameter<double>("DeltaIntersection") * CLHEP::mm;
57  m_stepMax = p.getParameter<double>("MaxStep") * CLHEP::cm;
58 
59  m_energyThreshold = p.getParameter<double>("EnergyThSimple") * CLHEP::GeV;
60 
61  // double
62  m_dChordSimple = p.getParameter<double>("DeltaChordSimple") * CLHEP::mm;
63  m_dOneStepSimple = p.getParameter<double>("DeltaOneStepSimple") * CLHEP::mm;
64  m_dIntersectionSimple = p.getParameter<double>("DeltaIntersectionSimple") * CLHEP::mm;
65  m_stepMaxSimple = p.getParameter<double>("MaxStepSimple") * CLHEP::cm;
66 
67  edm::LogVerbatim("SimG4CoreApplication")
68  << " New CMSFieldManager: LogicalVolume: <" << vol << ">\n"
69  << " Stepper: <" << stepper << ">\n"
70  << " Field type <" << type << ">\n"
71  << " Field const delta " << delta << " mm\n"
72  << " MaximumLoopCounts " << maxLC << "\n"
73  << " MinimumEpsilonStep " << minEpsStep << "\n"
74  << " MaximumEpsilonStep " << maxEpsStep << "\n"
75  << " MinStep " << minstep << " mm\n"
76  << " MaxStep " << m_stepMax / CLHEP::cm << " cm\n"
77  << " DeltaChord " << m_dChord << " mm\n"
78  << " DeltaOneStep " << m_dOneStep << " mm\n"
79  << " DeltaIntersection " << m_dIntersection << " mm\n"
80  << " EnergyThreshold " << m_energyThreshold << " MeV\n"
81  << " DeltaChordSimple " << m_dChordSimple << " mm\n"
82  << " DeltaOneStepSimple " << m_dOneStepSimple << " mm\n"
83  << " DeltaIntersectionSimple " << m_dIntersectionSimple << " mm\n"
84  << " MaxStepInVacuum " << m_stepMaxSimple / CLHEP::cm << " cm";
85 
86  // initialisation of chord finders
87  m_chordFinder = cf;
88  m_chordFinderMonopole = cfmon;
89 
90  // m_chordFinder->SetDeltaChord(dChord);
91  m_chordFinderMonopole->SetDeltaChord(m_dChord);
92 
93  // initialisation of field manager
94  theField.reset(field);
95  SetDetectorField(field);
96  SetMinimumEpsilonStep(minEpsStep);
97  SetMaximumEpsilonStep(maxEpsStep);
98 
99  // propagater in field
100  m_propagator = pf;
101  pf->SetMaxLoopCount(maxLC);
102  pf->SetMinimumEpsilonStep(minEpsStep);
103  pf->SetMaximumEpsilonStep(maxEpsStep);
104 
105  // initial initialisation the default chord finder is defined
106  SetMonopoleTracking(false);
107 
108  // define regions
109  std::vector<std::string> rnames = p.getParameter<std::vector<std::string>>("VacRegions");
110  if (!rnames.empty()) {
111  std::stringstream ss;
112  std::vector<G4Region *> *rs = G4RegionStore::GetInstance();
113  for (auto &regnam : rnames) {
114  for (auto &reg : *rs) {
115  if (regnam == reg->GetName()) {
116  m_regions.push_back(reg);
117  ss << " " << regnam;
118  }
119  }
120  }
121  edm::LogVerbatim("SimG4CoreApplication") << "Simple field integration in G4Regions:\n" << ss.str() << "\n";
122  }
123 }
124 
126  // run time parameters per track
127  if ((track->GetKineticEnergy() <= m_energyThreshold && track->GetParentID() > 0) || isInsideVacuum(track)) {
128  if (!m_cfVacuum) {
130  }
131 
132  } else if (m_cfVacuum) {
133  // restore defaults
135  }
136 }
137 
139  if (flag) {
141  if (m_cfVacuum) {
143  }
145  SetChordFinder(m_currChordFinder);
146  }
147  } else {
149  }
150  SetFieldChangesEnergy(flag);
151  m_currChordFinder->ResetStepEstimate();
152 }
153 
155  if (!m_regions.empty()) {
156  const G4Region *reg = track->GetVolume()->GetLogicalVolume()->GetRegion();
157  for (auto &areg : m_regions) {
158  if (reg == areg) {
159  return true;
160  }
161  }
162  }
163  return false;
164 }
165 
168  m_currChordFinder->SetDeltaChord(m_dChord);
169  SetChordFinder(m_currChordFinder);
170  SetDeltaOneStep(m_dOneStep);
171  SetDeltaIntersection(m_dIntersection);
172  m_propagator->SetLargestAcceptableStep(m_stepMax);
173  m_cfVacuum = false;
174 }
175 
177  m_currChordFinder->SetDeltaChord(m_dChordSimple);
178  SetDeltaOneStep(m_dOneStepSimple);
179  SetDeltaIntersection(m_dIntersectionSimple);
180  m_propagator->SetLargestAcceptableStep(m_stepMaxSimple);
181  m_cfVacuum = true;
182 }
CMSFieldManager::ConfigureForTrack
void ConfigureForTrack(const G4Track *) override
Definition: CMSFieldManager.cc:125
CMSFieldManager::m_dOneStepSimple
double m_dOneStepSimple
Definition: CMSFieldManager.h:68
CMSFieldManager::m_dOneStep
double m_dOneStep
Definition: CMSFieldManager.h:63
MessageLogger.h
funct::false
false
Definition: Factorize.h:34
CMSFieldManager::SetMonopoleTracking
void SetMonopoleTracking(G4bool)
Definition: CMSFieldManager.cc:138
CMSFieldManager::m_chordFinder
G4ChordFinder * m_chordFinder
Definition: CMSFieldManager.h:55
CMSFieldManager::m_dIntersectionSimple
double m_dIntersectionSimple
Definition: CMSFieldManager.h:69
CMSFieldManager::m_stepMax
double m_stepMax
Definition: CMSFieldManager.h:65
CMSFieldManager::m_stepMaxSimple
double m_stepMaxSimple
Definition: CMSFieldManager.h:70
CMSFieldManager.h
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:39
AlCaHLTBitMon_ParallelJobs.p
p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
CMSFieldManager::theField
std::unique_ptr< sim::Field > theField
Definition: CMSFieldManager.h:52
CMSFieldManager::isInsideVacuum
bool isInsideVacuum(const G4Track *)
Definition: CMSFieldManager.cc:154
contentValuesCheck.ss
ss
Definition: contentValuesCheck.py:33
CMSFieldManager::m_currChordFinder
G4ChordFinder * m_currChordFinder
Definition: CMSFieldManager.h:54
CMSFieldManager::m_chordFinderMonopole
G4ChordFinder * m_chordFinderMonopole
Definition: CMSFieldManager.h:56
sim::Field
Definition: Field.h:9
CMSFieldManager::m_dChord
double m_dChord
Definition: CMSFieldManager.h:62
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
CMSFieldManager::~CMSFieldManager
~CMSFieldManager() override
Definition: CMSFieldManager.cc:30
CMSFieldManager::m_propagator
G4PropagatorInField * m_propagator
Definition: CMSFieldManager.h:58
edm::ParameterSet
Definition: ParameterSet.h:36
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
edm::LogVerbatim
Definition: MessageLogger.h:297
CMSFieldManager::setChordFinderForVacuum
void setChordFinderForVacuum()
Definition: CMSFieldManager.cc:176
CMSFieldManager::setDefaultChordFinder
void setDefaultChordFinder()
Definition: CMSFieldManager.cc:166
type
type
Definition: HCALResponse.h:21
packedPFCandidateRefMixer_cfi.pf
pf
Definition: packedPFCandidateRefMixer_cfi.py:4
CMSFieldManager::m_dChordSimple
double m_dChordSimple
Definition: CMSFieldManager.h:67
HLT_2018_cff.track
track
Definition: HLT_2018_cff.py:10352
CMSFieldManager::CMSFieldManager
CMSFieldManager()
Definition: CMSFieldManager.cc:13
CMSFieldManager::m_regions
std::vector< const G4Region * > m_regions
Definition: CMSFieldManager.h:60
CMSFieldManager::m_energyThreshold
double m_energyThreshold
Definition: CMSFieldManager.h:66
CMSFieldManager::m_dIntersection
double m_dIntersection
Definition: CMSFieldManager.h:64
RemoveAddSevLevel.flag
flag
Definition: RemoveAddSevLevel.py:116
CMSFieldManager::m_cfVacuum
bool m_cfVacuum
Definition: CMSFieldManager.h:72