CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
FieldBuilder.cc
Go to the documentation of this file.
1 // #include "DetectorDescription/Core/interface/DDLogicalPart.h"
2 
5 
8 
12 #include "SimG4Core/MagneticField/interface/G4MonopoleEquation.hh"
14 
15 /*
16 #include "DetectorDescription/Parser/interface/DDLConfiguration.h"
17 #include "DetectorDescription/Base/interface/DDException.h"
18 #include "DetectorDescription/Algorithm/src/AlgoInit.h"
19 */
20 
21 #include "G4Mag_UsualEqRhs.hh"
22 #include "G4ClassicalRK4.hh"
23 #include "G4PropagatorInField.hh"
24 #include "G4FieldManager.hh"
25 #include "G4TransportationManager.hh"
26 #include "G4ChordFinder.hh"
27 #include "G4UniformMagField.hh"
28 
30 
31 #include "G4LogicalVolumeStore.hh"
32 
33 using namespace sim;
34 
36  const edm::ParameterSet & p)
37  : theField( new Field(f,p.getParameter<double>("delta"))),
38  theFieldEquation(new G4Mag_UsualEqRhs(theField.get())),
39  theTopVolume(0), fChordFinder(0), fChordFinderMonopole(0),
40  fieldValue(0.), minStep(0.), dChord(0.), dOneStep(0.),
41  dIntersection(0.), dIntersectionAndOneStep(0.),
42  maxLoopCount(0), minEpsilonStep(0.), maxEpsilonStep(0.),
43  thePSet(p) {
44 
45  theField->fieldEquation(theFieldEquation);
46 }
47 
48 
49 void FieldBuilder::build( G4FieldManager* fM, G4PropagatorInField* fP) {
50 
51  edm::ParameterSet thePSetForGMFM =
52  thePSet.getParameter<edm::ParameterSet>("ConfGlobalMFM");
53 
54  std::string volName = thePSetForGMFM.getParameter< std::string >("Volume");
55 
56  edm::ParameterSet volPSet =
57  thePSetForGMFM.getParameter< edm::ParameterSet >( volName );
58 
59  configureForVolume( volName, volPSet, fM, fP );
60 
61  // configure( "MagneticFieldType", fM, fP ) ;
62 
63  if ( thePSet.getParameter<bool>("UseLocalMagFieldManager") ) {
64 
65  edm::ParameterSet defpset ;
66  edm::ParameterSet thePSetForLMFM =
67  thePSet.getUntrackedParameter<edm::ParameterSet>("ConfLocalMFM", defpset);
68  //
69  // Patology !!! LocalFM requested but configuration not given !
70  // In principal, need to throw an exception
71  //
72  if ( thePSetForLMFM == defpset ) {
73  std::cout << " Patology ! Local Mag.Field Manager requested but config not given !\n";
74  return ;
75  }
76 
77  std::vector<std::string> ListOfVolumes =
78  thePSetForLMFM.getParameter< std::vector<std::string> >("ListOfVolumes");
79 
80  // creating Local Mag.Field Manager
81  for (unsigned int i = 0; i < ListOfVolumes.size(); ++ i ) {
82  volPSet = thePSetForLMFM.getParameter< edm::ParameterSet >(ListOfVolumes[i]);
83  G4FieldManager* fAltM = new G4FieldManager() ;
84  configureForVolume( ListOfVolumes[i], volPSet, fAltM ) ;
85  //configureLocalFM( ListOfVolumes[i], fAltM ) ;
86  LocalFieldManager* fLM = new LocalFieldManager( theField.get(), fM, fAltM ) ;
87  fLM->SetVerbosity(thePSet.getUntrackedParameter<bool>("Verbosity",false));
88  theTopVolume->SetFieldManager( fLM, true ) ;
89  }
90  }
91  return ;
92 
93 }
94 
96  edm::ParameterSet& volPSet,
97  G4FieldManager * fM,
98  G4PropagatorInField * fP ) {
99 
100  G4LogicalVolumeStore* theStore = G4LogicalVolumeStore::GetInstance();
101  for (unsigned int i=0; i<(*theStore).size(); ++i ) {
102  std::string curVolName = ((*theStore)[i])->GetName();
103  if ( curVolName == volName ) {
104  theTopVolume = (*theStore)[i] ;
105  }
106  }
107 
108  fieldType = volPSet.getParameter<std::string>("Type") ;
109  stepper = volPSet.getParameter<std::string>("Stepper") ;
110  edm::ParameterSet stpPSet =
112  minStep = stpPSet.getParameter<double>("MinStep") ;
113  dChord = stpPSet.getParameter<double>("DeltaChord") ;
114  dOneStep = stpPSet.getParameter<double>("DeltaOneStep") ;
115  dIntersection = stpPSet.getParameter<double>("DeltaIntersection") ;
116  dIntersectionAndOneStep = stpPSet.getUntrackedParameter<double>("DeltaIntersectionAndOneStep",-1.);
117  maxLoopCount = stpPSet.getUntrackedParameter<double>("MaximumLoopCounts",1000);
118  minEpsilonStep = stpPSet.getUntrackedParameter<double>("MinimumEpsilonStep",0.00001);
119  maxEpsilonStep = stpPSet.getUntrackedParameter<double>("MaximumEpsilonStep",0.01);
120 
121  if (fM!=0) configureFieldManager(fM);
122  if (fP!=0) configurePropagatorInField(fP);
123 
124  return;
125 
126 }
127 
128 G4LogicalVolume * FieldBuilder::fieldTopVolume() { return theTopVolume; }
129 
130 void FieldBuilder::setStepperAndChordFinder(G4FieldManager * fM, int val) {
131 
132  if (fM != 0) {
133  if (val == 0) {
134  if (fChordFinder != 0) fM->SetChordFinder(fChordFinder);
135  } else {
136  fChordFinder = fM->GetChordFinder();
137  if (fChordFinderMonopole != 0) fM->SetChordFinder(fChordFinderMonopole);
138  }
139  }
140 }
141 
142 
143 void FieldBuilder::configureFieldManager(G4FieldManager * fM) {
144 
145  if (fM!=0) {
146  fM->SetDetectorField(theField.get());
147  FieldStepper * theStepper = new FieldStepper(theField->fieldEquation());
148  theStepper->select(stepper);
149  G4ChordFinder * CF = new G4ChordFinder(theField.get(),minStep,theStepper);
150  CF->SetDeltaChord(dChord);
151  fM->SetChordFinder(CF);
152  fM->SetDeltaOneStep(dOneStep);
153  fM->SetDeltaIntersection(dIntersection);
154  if (dIntersectionAndOneStep != -1.)
155  fM->SetAccuraciesWithDeltaOneStep(dIntersectionAndOneStep);
156  }
157  if (fChordFinderMonopole == 0) {
158  G4MonopoleEquation* fMonopoleEquation = new G4MonopoleEquation(theField.get());
159  G4MagIntegratorStepper* theStepper = new G4ClassicalRK4(fMonopoleEquation,8);
160  fChordFinderMonopole = new G4ChordFinder(theField.get(),minStep,theStepper);
161  fChordFinderMonopole->SetDeltaChord(dChord);
162  }
163 }
164 
165 void FieldBuilder::configurePropagatorInField(G4PropagatorInField * fP) {
166  if (fP==0) return;
167  fP->SetMaxLoopCount(int(maxLoopCount));
168  fP->SetMinimumEpsilonStep(minEpsilonStep);
169  fP->SetMaximumEpsilonStep(maxEpsilonStep);
170  fP->SetVerboseLevel(0);
171  return;
172 }
173 
T getParameter(std::string const &) const
G4ChordFinder * fChordFinder
Definition: FieldBuilder.h:50
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
std::string fieldType
Definition: FieldBuilder.h:53
void build(G4FieldManager *fM=0, G4PropagatorInField *fP=0)
Definition: FieldBuilder.cc:49
G4ChordFinder * fChordFinderMonopole
Definition: FieldBuilder.h:50
Definition: sim.h:19
void configureForVolume(const std::string &volName, edm::ParameterSet &volPSet, G4FieldManager *fM=0, G4PropagatorInField *fP=0)
Definition: FieldBuilder.cc:95
G4LogicalVolume * theTopVolume
Definition: FieldBuilder.h:49
void SetVerbosity(bool flag)
double f[11][100]
void setStepperAndChordFinder(G4FieldManager *fM, int val)
double dIntersectionAndOneStep
Definition: FieldBuilder.h:60
FieldBuilder(const MagneticField *, const edm::ParameterSet &)
Definition: FieldBuilder.cc:35
G4MagIntegratorStepper * select(const std::string &s)
Definition: FieldStepper.cc:32
void configureFieldManager(G4FieldManager *fM)
std::auto_ptr< Field > theField
Definition: FieldBuilder.h:47
G4LogicalVolume * fieldTopVolume()
G4Mag_UsualEqRhs * theFieldEquation
Definition: FieldBuilder.h:48
tuple cout
Definition: gather_cfg.py:121
edm::ParameterSet thePSet
Definition: FieldBuilder.h:64
std::string stepper
Definition: FieldBuilder.h:55
void configurePropagatorInField(G4PropagatorInField *fP)
T get(const Candidate &c)
Definition: component.h:56