CMS 3D CMS Logo

CMSFieldManager.cc
Go to the documentation of this file.
4 
5 #include "G4ChordFinder.hh"
6 #include "G4Track.hh"
7 #include "CLHEP/Units/GlobalSystemOfUnits.h"
8 
10  : G4FieldManager(), currChordFinder(nullptr), chordFinder(nullptr),
11  chordFinderMonopole(nullptr), dChord(0.001), dOneStep(0.001),
12  dIntersection(0.0001), energyThreshold(0.0), dChordSimple(0.1),
13  dOneStepSimple(0.1), dIntersectionSimple(0.01)
14 {}
15 
17 {}
18 
20  G4ChordFinder* cf, G4ChordFinder* cfmon,
21  const std::string& vol, const std::string& type,
22  const std::string& stepper,
23  double delta, double minstep)
24 {
25  dChord = p.getParameter<double>("DeltaChord")*CLHEP::mm;
26  dOneStep = p.getParameter<double>("DeltaOneStep")*CLHEP::mm;
27  dIntersection = p.getParameter<double>("DeltaIntersection")*CLHEP::mm;
28  energyThreshold = p.getParameter<double>("EnergyThSimple")*CLHEP::GeV;
29  dChordSimple = p.getParameter<double>("DeltaChordSimple")*CLHEP::mm;
30  dOneStepSimple = p.getParameter<double>("DeltaOneStepSimple")*CLHEP::mm;
31  dIntersectionSimple = p.getParameter<double>("DeltaIntersectionSimple")*CLHEP::mm;
32  int maxLC = (int)p.getUntrackedParameter<double>("MaximumLoopCounts",1000);
33  double minEpsStep =
34  p.getUntrackedParameter<double>("MinimumEpsilonStep",0.00001)*CLHEP::mm;
35  double maxEpsStep =
36  p.getUntrackedParameter<double>("MaximumEpsilonStep",0.01)*CLHEP::mm;
37  edm::LogInfo("SimG4CoreApplication")
38  << " New CMSFieldManager: LogicalVolume: <" << vol << ">\n"
39  << " Stepper: <" << stepper << ">\n"
40  << " Field type <" << type<< ">\n"
41  << " Field const delta(mm) " << delta << "\n"
42  << " MinStep(mm) " << minstep<< "\n"
43  << " DeltaChord(mm) " << dChord<< "\n"
44  << " DeltaOneStep(mm) " << dOneStep<< "\n"
45  << " DeltaIntersection(mm) " << dIntersection<< "\n"
46  << " EnergyThreshold(GeV) " << energyThreshold<< "\n"
47  << " DeltaChordSimple(mm) " << dChord<< "\n"
48  << " DeltaOneStepSimple(mm) " << dOneStep<< "\n"
49  << " DeltaIntersectionSimple(mm) " << dIntersection<< "\n"
50  << " MaximumLoopCounts " << maxLC<< "\n"
51  << " MinimumEpsilonStep " << minEpsStep<< "\n"
52  << " MaximumEpsilonStep " << maxEpsStep;
53 
54  // initialisation of chord finders
55  chordFinder = cf;
56  chordFinderMonopole = cfmon;
57  chordFinder->SetDeltaChord(dChord);
58  chordFinderMonopole->SetDeltaChord(dChord);
59 
60  // initialisation of field manager
61  theField.reset(field);
62  SetDetectorField(field);
63  SetDeltaOneStep(dOneStep);
64  SetDeltaIntersection(dIntersection);
65  SetMinimumEpsilonStep(minEpsStep);
66  SetMaximumEpsilonStep(maxEpsStep);
67 
68  SetMonopoleTracking(false);
69 }
70 
72 {
73  // run time parameters per track
74  if(track->GetKineticEnergy() <= energyThreshold && track->GetParentID() > 0) {
75  chordFinder->SetDeltaChord(dChordSimple);
76  SetDeltaOneStep(dOneStepSimple);
77  SetDeltaIntersection(dIntersectionSimple);
78  } else {
79  chordFinder->SetDeltaChord(dChord);
80  SetDeltaOneStep(dOneStep);
81  SetDeltaIntersection(dIntersection);
82  }
83 }
84 
86 {
87  if(flag) {
89  SetFieldChangesEnergy(true);
90  } else {
92  SetFieldChangesEnergy(false);
93  }
94  SetChordFinder(currChordFinder);
95 }
dbl * delta
Definition: mlp_gen.cc:36
type
Definition: HCALResponse.h:21
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
const double GeV
Definition: MathUtil.h:16
void InitialiseForVolume(const edm::ParameterSet &, sim::Field *, G4ChordFinder *, G4ChordFinder *, const std::string &vol, const std::string &type, const std::string &stepper, double delta, double minstep)
G4ChordFinder * chordFinder
G4ChordFinder * chordFinderMonopole
#define nullptr
double dIntersectionSimple
G4ChordFinder * currChordFinder
~CMSFieldManager() override
void ConfigureForTrack(const G4Track *) override
std::unique_ptr< sim::Field > theField
void SetMonopoleTracking(G4bool)