CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
CMSFieldManager Class Reference

#include <CMSFieldManager.h>

Inheritance diagram for CMSFieldManager:

Public Member Functions

 CMSFieldManager ()
 
 CMSFieldManager (const CMSFieldManager &)=delete
 
void ConfigureForTrack (const G4Track *) override
 
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 *)
 
CMSFieldManageroperator= (const CMSFieldManager &)=delete
 
void setMonopoleTracking (G4bool)
 
 ~CMSFieldManager () override
 

Private Member Functions

bool isInsideTracker (const G4Track *)
 
bool isInsideVacuum (const G4Track *)
 
void setChordFinderForTracker ()
 
void setChordFinderForVacuum ()
 
void setDefaultChordFinder ()
 

Private Attributes

bool m_cfTracker
 
bool m_cfVacuum
 
G4ChordFinder * m_chordFinder
 
G4ChordFinder * m_chordFinderMonopole
 
G4ChordFinder * m_currChordFinder
 
double m_dChord
 
double m_dChordSimple
 
double m_dChordTracker
 
double m_dIntersection
 
double m_dIntersectionSimple
 
double m_dInterTracker
 
double m_dOneStep
 
double m_dOneStepSimple
 
double m_dOneStepTracker
 
double m_energyThreshold
 
double m_energyThTracker
 
G4PropagatorInField * m_propagator
 
std::vector< const G4Region * > m_regions
 
double m_Rmax2
 
double m_stepMax
 
double m_stepMaxSimple
 
double m_Zmax
 
std::unique_ptr< sim::FieldtheField
 

Detailed Description

Definition at line 24 of file CMSFieldManager.h.

Constructor & Destructor Documentation

◆ CMSFieldManager() [1/2]

CMSFieldManager::CMSFieldManager ( )
explicit

Definition at line 13 of file CMSFieldManager.cc.

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),
33  m_stepMaxSimple(1000.),
34  m_cfTracker(false),
35  m_cfVacuum(false) {}
G4ChordFinder * m_chordFinder
double m_dIntersectionSimple
double m_energyThreshold
double m_dOneStepTracker
G4PropagatorInField * m_propagator
G4ChordFinder * m_currChordFinder
G4ChordFinder * m_chordFinderMonopole
double m_energyThTracker

◆ ~CMSFieldManager()

CMSFieldManager::~CMSFieldManager ( )
override

Definition at line 37 of file CMSFieldManager.cc.

References m_chordFinder, m_chordFinderMonopole, and m_currChordFinder.

37  {
39  delete m_chordFinder;
40  }
42  delete m_chordFinderMonopole;
43  }
44 }
G4ChordFinder * m_chordFinder
G4ChordFinder * m_currChordFinder
G4ChordFinder * m_chordFinderMonopole

◆ CMSFieldManager() [2/2]

CMSFieldManager::CMSFieldManager ( const CMSFieldManager )
delete

Member Function Documentation

◆ ConfigureForTrack()

void CMSFieldManager::ConfigureForTrack ( const G4Track *  track)
override

Definition at line 140 of file CMSFieldManager.cc.

References isInsideTracker(), isInsideVacuum(), m_cfTracker, m_cfVacuum, m_energyThreshold, m_energyThTracker, setChordFinderForTracker(), setChordFinderForVacuum(), setDefaultChordFinder(), and HLT_2024v13_cff::track.

Referenced by MonopoleTransportation::AlongStepGetPhysicalInteractionLength().

140  {
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 }
void setDefaultChordFinder()
double m_energyThreshold
void setChordFinderForTracker()
bool isInsideVacuum(const G4Track *)
void setChordFinderForVacuum()
double m_energyThTracker
bool isInsideTracker(const G4Track *)

◆ InitialiseForVolume()

void CMSFieldManager::InitialiseForVolume ( const edm::ParameterSet p,
sim::Field field,
G4ChordFinder *  cfDefault,
G4ChordFinder *  cfMonopole,
const std::string &  vol,
const std::string &  fieldType,
const std::string &  stepperName,
double  delta,
G4PropagatorInField *  pf 
)

Definition at line 46 of file CMSFieldManager.cc.

References dumpMFGeometry_cfg::delta, createfilelist::int, m_chordFinder, m_chordFinderMonopole, m_dChord, m_dChordSimple, m_dChordTracker, m_dIntersection, m_dIntersectionSimple, m_dInterTracker, m_dOneStep, m_dOneStepSimple, m_dOneStepTracker, m_energyThreshold, m_energyThTracker, m_propagator, m_regions, m_Rmax2, m_stepMax, m_stepMaxSimple, m_Zmax, AlCaHLTBitMon_ParallelJobs::p, packedPFCandidateRefMixer_cfi::pf, setMonopoleTracking(), contentValuesCheck::ss, and theField.

Referenced by sim::FieldBuilder::configureForVolume().

54  {
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 }
G4ChordFinder * m_chordFinder
Log< level::Info, true > LogVerbatim
double m_dIntersectionSimple
double m_energyThreshold
double m_dOneStepTracker
G4PropagatorInField * m_propagator
G4ChordFinder * m_chordFinderMonopole
std::unique_ptr< sim::Field > theField
std::vector< const G4Region * > m_regions
void setMonopoleTracking(G4bool)
double m_energyThTracker

◆ isInsideTracker()

bool CMSFieldManager::isInsideTracker ( const G4Track *  track)
private

Definition at line 186 of file CMSFieldManager.cc.

References funct::abs(), m_Rmax2, m_Zmax, HLT_2024v13_cff::track, x, and y.

Referenced by ConfigureForTrack().

186  {
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 }
Abs< T >::type abs(const T &t)
Definition: Abs.h:22

◆ isInsideVacuum()

bool CMSFieldManager::isInsideVacuum ( const G4Track *  track)
private

Definition at line 174 of file CMSFieldManager.cc.

References m_regions, and HLT_2024v13_cff::track.

Referenced by ConfigureForTrack().

174  {
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 }
std::vector< const G4Region * > m_regions

◆ operator=()

CMSFieldManager& CMSFieldManager::operator= ( const CMSFieldManager )
delete

◆ setChordFinderForTracker()

void CMSFieldManager::setChordFinderForTracker ( )
private

Definition at line 205 of file CMSFieldManager.cc.

References m_cfTracker, m_cfVacuum, m_chordFinder, m_currChordFinder, m_dChordTracker, m_dInterTracker, m_dOneStepTracker, m_propagator, and m_stepMax.

Referenced by ConfigureForTrack().

205  {
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 }
G4ChordFinder * m_chordFinder
double m_dOneStepTracker
G4PropagatorInField * m_propagator
G4ChordFinder * m_currChordFinder

◆ setChordFinderForVacuum()

void CMSFieldManager::setChordFinderForVacuum ( )
private

Definition at line 218 of file CMSFieldManager.cc.

References m_cfTracker, m_cfVacuum, m_chordFinder, m_currChordFinder, m_dChordSimple, m_dIntersectionSimple, m_dOneStepSimple, m_propagator, and m_stepMaxSimple.

Referenced by ConfigureForTrack().

218  {
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 }
G4ChordFinder * m_chordFinder
double m_dIntersectionSimple
G4PropagatorInField * m_propagator
G4ChordFinder * m_currChordFinder

◆ setDefaultChordFinder()

void CMSFieldManager::setDefaultChordFinder ( )
private

Definition at line 193 of file CMSFieldManager.cc.

References m_cfTracker, m_cfVacuum, m_chordFinder, m_currChordFinder, m_dChord, m_dIntersection, m_dOneStep, m_propagator, and m_stepMax.

Referenced by ConfigureForTrack(), and setMonopoleTracking().

193  {
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 }
G4ChordFinder * m_chordFinder
G4PropagatorInField * m_propagator
G4ChordFinder * m_currChordFinder

◆ setMonopoleTracking()

void CMSFieldManager::setMonopoleTracking ( G4bool  flag)

Definition at line 158 of file CMSFieldManager.cc.

References RemoveAddSevLevel::flag, m_cfTracker, m_cfVacuum, m_chordFinderMonopole, m_currChordFinder, and setDefaultChordFinder().

Referenced by MonopoleTransportation::AlongStepGetPhysicalInteractionLength(), InitialiseForVolume(), and MonopoleTransportation::PostStepDoIt().

158  {
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 }
void setDefaultChordFinder()
G4ChordFinder * m_currChordFinder
G4ChordFinder * m_chordFinderMonopole

Member Data Documentation

◆ m_cfTracker

bool CMSFieldManager::m_cfTracker
private

◆ m_cfVacuum

bool CMSFieldManager::m_cfVacuum
private

◆ m_chordFinder

G4ChordFinder* CMSFieldManager::m_chordFinder
private

◆ m_chordFinderMonopole

G4ChordFinder* CMSFieldManager::m_chordFinderMonopole
private

Definition at line 58 of file CMSFieldManager.h.

Referenced by InitialiseForVolume(), setMonopoleTracking(), and ~CMSFieldManager().

◆ m_currChordFinder

G4ChordFinder* CMSFieldManager::m_currChordFinder
private

◆ m_dChord

double CMSFieldManager::m_dChord
private

Definition at line 64 of file CMSFieldManager.h.

Referenced by InitialiseForVolume(), and setDefaultChordFinder().

◆ m_dChordSimple

double CMSFieldManager::m_dChordSimple
private

Definition at line 75 of file CMSFieldManager.h.

Referenced by InitialiseForVolume(), and setChordFinderForVacuum().

◆ m_dChordTracker

double CMSFieldManager::m_dChordTracker
private

Definition at line 65 of file CMSFieldManager.h.

Referenced by InitialiseForVolume(), and setChordFinderForTracker().

◆ m_dIntersection

double CMSFieldManager::m_dIntersection
private

Definition at line 68 of file CMSFieldManager.h.

Referenced by InitialiseForVolume(), and setDefaultChordFinder().

◆ m_dIntersectionSimple

double CMSFieldManager::m_dIntersectionSimple
private

Definition at line 77 of file CMSFieldManager.h.

Referenced by InitialiseForVolume(), and setChordFinderForVacuum().

◆ m_dInterTracker

double CMSFieldManager::m_dInterTracker
private

Definition at line 69 of file CMSFieldManager.h.

Referenced by InitialiseForVolume(), and setChordFinderForTracker().

◆ m_dOneStep

double CMSFieldManager::m_dOneStep
private

Definition at line 66 of file CMSFieldManager.h.

Referenced by InitialiseForVolume(), and setDefaultChordFinder().

◆ m_dOneStepSimple

double CMSFieldManager::m_dOneStepSimple
private

Definition at line 76 of file CMSFieldManager.h.

Referenced by InitialiseForVolume(), and setChordFinderForVacuum().

◆ m_dOneStepTracker

double CMSFieldManager::m_dOneStepTracker
private

Definition at line 67 of file CMSFieldManager.h.

Referenced by InitialiseForVolume(), and setChordFinderForTracker().

◆ m_energyThreshold

double CMSFieldManager::m_energyThreshold
private

Definition at line 74 of file CMSFieldManager.h.

Referenced by ConfigureForTrack(), and InitialiseForVolume().

◆ m_energyThTracker

double CMSFieldManager::m_energyThTracker
private

Definition at line 73 of file CMSFieldManager.h.

Referenced by ConfigureForTrack(), and InitialiseForVolume().

◆ m_propagator

G4PropagatorInField* CMSFieldManager::m_propagator
private

◆ m_regions

std::vector<const G4Region *> CMSFieldManager::m_regions
private

Definition at line 62 of file CMSFieldManager.h.

Referenced by InitialiseForVolume(), and isInsideVacuum().

◆ m_Rmax2

double CMSFieldManager::m_Rmax2
private

Definition at line 70 of file CMSFieldManager.h.

Referenced by InitialiseForVolume(), and isInsideTracker().

◆ m_stepMax

double CMSFieldManager::m_stepMax
private

◆ m_stepMaxSimple

double CMSFieldManager::m_stepMaxSimple
private

Definition at line 78 of file CMSFieldManager.h.

Referenced by InitialiseForVolume(), and setChordFinderForVacuum().

◆ m_Zmax

double CMSFieldManager::m_Zmax
private

Definition at line 71 of file CMSFieldManager.h.

Referenced by InitialiseForVolume(), and isInsideTracker().

◆ theField

std::unique_ptr<sim::Field> CMSFieldManager::theField
private

Definition at line 54 of file CMSFieldManager.h.

Referenced by InitialiseForVolume().