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