CMS 3D CMS Logo

FieldBuilder.cc
Go to the documentation of this file.
3 
5 
11 
12 #include "CLHEP/Units/GlobalSystemOfUnits.h"
13 #include "G4ChordFinder.hh"
14 #include "G4ClassicalRK4.hh"
15 #include "G4FieldManager.hh"
16 #include "G4LogicalVolumeStore.hh"
17 #include "G4Mag_UsualEqRhs.hh"
18 #include "G4TMagFieldEquation.hh"
19 #include "G4PropagatorInField.hh"
20 
21 using namespace sim;
22 
23 FieldBuilder::FieldBuilder(const MagneticField *f, const edm::ParameterSet &p) : theTopVolume(nullptr), thePSet(p) {
24  theDelta = p.getParameter<double>("delta") * CLHEP::mm;
25  theField = new Field(f, theDelta);
26  theFieldEquation = nullptr;
27 }
28 
30 
31 void FieldBuilder::build(CMSFieldManager *fM, G4PropagatorInField *fP) {
32  edm::ParameterSet thePSetForGMFM = thePSet.getParameter<edm::ParameterSet>("ConfGlobalMFM");
33  std::string volName = thePSetForGMFM.getParameter<std::string>("Volume");
34  edm::ParameterSet volPSet = thePSetForGMFM.getParameter<edm::ParameterSet>(volName);
35 
36  configureForVolume(volName, volPSet, fM, fP);
37 
38  edm::LogVerbatim("SimG4CoreMagneticField") << " FieldBuilder::build: Global magnetic field is used";
39 }
40 
42  edm::ParameterSet &volPSet,
43  CMSFieldManager *fM,
44  G4PropagatorInField *fP) {
45  G4LogicalVolumeStore *theStore = G4LogicalVolumeStore::GetInstance();
46  for (auto vol : *theStore) {
47  if ((std::string)vol->GetName() == volName) {
48  theTopVolume = vol;
49  break;
50  }
51  }
52 
53  std::string fieldType = volPSet.getParameter<std::string>("Type");
54  std::string stepper = volPSet.getParameter<std::string>("Stepper");
55 
56  edm::ParameterSet stpPSet = volPSet.getParameter<edm::ParameterSet>("StepperParam");
57  double minStep = stpPSet.getParameter<double>("MinStep") * CLHEP::mm;
58 
59  if (stepper == "G4TDormandPrince45") {
60  theFieldEquation = new G4TMagFieldEquation<Field>(theField);
61  } else {
62  theFieldEquation = new G4Mag_UsualEqRhs(theField);
63  }
64 
65  FieldStepper *dStepper = new FieldStepper(theFieldEquation, theDelta, stepper);
66  G4ChordFinder *cf = new G4ChordFinder(theField, minStep, dStepper);
67 
68  MonopoleEquation *monopoleEquation = new MonopoleEquation(theField);
69  G4MagIntegratorStepper *mStepper = new G4ClassicalRK4(monopoleEquation, 8);
70  G4ChordFinder *cfmon = new G4ChordFinder(theField, minStep, mStepper);
71 
72  fM->InitialiseForVolume(stpPSet, theField, cf, cfmon, volName, fieldType, stepper, theDelta, fP);
73 }
Log< level::Info, true > LogVerbatim
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
G4LogicalVolume * theTopVolume
Definition: FieldBuilder.h:31
double f[11][100]
FieldBuilder(const MagneticField *, const edm::ParameterSet &)
Definition: FieldBuilder.cc:23
void build(CMSFieldManager *fM, G4PropagatorInField *fP)
Definition: FieldBuilder.cc:31
G4Mag_UsualEqRhs * theFieldEquation
Definition: FieldBuilder.h:30
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 *)
edm::ParameterSet thePSet
Definition: FieldBuilder.h:32
void configureForVolume(const std::string &volName, edm::ParameterSet &volPSet, CMSFieldManager *fM, G4PropagatorInField *fP)
Definition: FieldBuilder.cc:41