CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_13_patch3/src/SimG4Core/MagneticField/src/FieldBuilder.cc

Go to the documentation of this file.
00001 // #include "DetectorDescription/Core/interface/DDLogicalPart.h"
00002 
00003 #include "MagneticField/Engine/interface/MagneticField.h"
00004 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
00005 
00006 #include "FWCore/Framework/interface/ESHandle.h"
00007 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00008 
00009 #include "SimG4Core/MagneticField/interface/FieldBuilder.h"
00010 #include "SimG4Core/MagneticField/interface/Field.h"
00011 #include "SimG4Core/MagneticField/interface/FieldStepper.h"
00012 #include "SimG4Core/MagneticField/interface/G4MonopoleEquation.hh"
00013 #include "SimG4Core/Notification/interface/SimG4Exception.h"
00014 
00015 /*
00016 #include "DetectorDescription/Parser/interface/DDLConfiguration.h"
00017 #include "DetectorDescription/Base/interface/DDException.h"
00018 #include "DetectorDescription/Algorithm/src/AlgoInit.h"
00019 */
00020 
00021 #include "G4Mag_UsualEqRhs.hh"
00022 #include "G4ClassicalRK4.hh"
00023 #include "G4PropagatorInField.hh"
00024 #include "G4FieldManager.hh"
00025 #include "G4TransportationManager.hh"
00026 #include "G4ChordFinder.hh"
00027 #include "G4UniformMagField.hh"
00028 
00029 #include "SimG4Core/MagneticField/interface/LocalFieldManager.h"
00030 
00031 #include "G4LogicalVolumeStore.hh"
00032 
00033 using namespace sim;
00034 
00035 FieldBuilder::FieldBuilder(const MagneticField * f, 
00036                            const edm::ParameterSet & p) 
00037    : theField( new Field(f,p.getParameter<double>("delta"))), 
00038      theFieldEquation(new G4Mag_UsualEqRhs(theField.get())),
00039      theTopVolume(0), fChordFinder(0), fChordFinderMonopole(0),
00040      fieldValue(0.), minStep(0.), dChord(0.), dOneStep(0.),
00041      dIntersection(0.), dIntersectionAndOneStep(0.), 
00042      maxLoopCount(0), minEpsilonStep(0.), maxEpsilonStep(0.), 
00043      thePSet(p) {
00044 
00045   theField->fieldEquation(theFieldEquation);
00046 }
00047 
00048 
00049 void FieldBuilder::build( G4FieldManager* fM, G4PropagatorInField* fP) {
00050     
00051   edm::ParameterSet thePSetForGMFM =
00052     thePSet.getParameter<edm::ParameterSet>("ConfGlobalMFM");
00053 
00054   std::string volName = thePSetForGMFM.getParameter< std::string >("Volume");
00055   
00056   edm::ParameterSet volPSet =
00057     thePSetForGMFM.getParameter< edm::ParameterSet >( volName );
00058     
00059   configureForVolume( volName, volPSet, fM, fP );
00060     
00061   // configure( "MagneticFieldType", fM, fP ) ;
00062 
00063   if ( thePSet.getParameter<bool>("UseLocalMagFieldManager") )  {
00064 
00065     edm::ParameterSet defpset ;
00066     edm::ParameterSet thePSetForLMFM = 
00067       thePSet.getUntrackedParameter<edm::ParameterSet>("ConfLocalMFM", defpset);
00068     //
00069     // Patology !!! LocalFM requested but configuration not given ! 
00070     // In principal, need to throw an exception
00071     //
00072     if ( thePSetForLMFM == defpset )  {
00073       std::cout << " Patology ! Local Mag.Field Manager requested but config not given !\n";
00074       return ;
00075     }
00076        
00077     std::vector<std::string> ListOfVolumes = 
00078       thePSetForLMFM.getParameter< std::vector<std::string> >("ListOfVolumes");
00079           
00080     // creating Local Mag.Field Manager
00081     for (unsigned int i = 0; i < ListOfVolumes.size(); ++ i )   {
00082       volPSet = thePSetForLMFM.getParameter< edm::ParameterSet >(ListOfVolumes[i]);
00083       G4FieldManager* fAltM = new G4FieldManager() ;
00084       configureForVolume( ListOfVolumes[i], volPSet, fAltM ) ;
00085       //configureLocalFM( ListOfVolumes[i], fAltM ) ;
00086       LocalFieldManager* fLM = new LocalFieldManager( theField.get(), fM, fAltM ) ;
00087       fLM->SetVerbosity(thePSet.getUntrackedParameter<bool>("Verbosity",false));
00088       theTopVolume->SetFieldManager( fLM, true ) ;
00089     }
00090   }
00091   return ;
00092  
00093 }
00094 
00095 void FieldBuilder::configureForVolume( const std::string& volName,
00096                                        edm::ParameterSet& volPSet,
00097                                        G4FieldManager * fM,
00098                                        G4PropagatorInField * fP ) {
00099 
00100   G4LogicalVolumeStore* theStore = G4LogicalVolumeStore::GetInstance();
00101   for (unsigned int i=0; i<(*theStore).size(); ++i ) {
00102     std::string curVolName = ((*theStore)[i])->GetName();
00103     if ( curVolName == volName ) {
00104       theTopVolume = (*theStore)[i] ;
00105     }
00106   }
00107 
00108   fieldType     = volPSet.getParameter<std::string>("Type") ;
00109   stepper       = volPSet.getParameter<std::string>("Stepper") ;
00110   edm::ParameterSet stpPSet = 
00111     volPSet.getParameter<edm::ParameterSet>(stepper) ;
00112   minStep       = stpPSet.getParameter<double>("MinStep") ;
00113   dChord        = stpPSet.getParameter<double>("DeltaChord") ;
00114   dOneStep      = stpPSet.getParameter<double>("DeltaOneStep") ;
00115   dIntersection = stpPSet.getParameter<double>("DeltaIntersection") ;
00116   dIntersectionAndOneStep = stpPSet.getUntrackedParameter<double>("DeltaIntersectionAndOneStep",-1.);
00117   maxLoopCount = stpPSet.getUntrackedParameter<double>("MaximumLoopCounts",1000);
00118   minEpsilonStep = stpPSet.getUntrackedParameter<double>("MinimumEpsilonStep",0.00001);
00119   maxEpsilonStep = stpPSet.getUntrackedParameter<double>("MaximumEpsilonStep",0.01);
00120    
00121   if (fM!=0) configureFieldManager(fM);
00122   if (fP!=0) configurePropagatorInField(fP);    
00123 
00124   return;
00125 
00126 }
00127 
00128 G4LogicalVolume * FieldBuilder::fieldTopVolume() { return theTopVolume; }
00129 
00130 void FieldBuilder::setStepperAndChordFinder(G4FieldManager * fM, int val) {
00131 
00132   if (fM != 0) {
00133     if (val == 0) {
00134       if (fChordFinder != 0) fM->SetChordFinder(fChordFinder);
00135     } else {
00136       fChordFinder = fM->GetChordFinder();
00137       if (fChordFinderMonopole != 0) fM->SetChordFinder(fChordFinderMonopole);
00138     }
00139   }
00140 }
00141 
00142  
00143 void FieldBuilder::configureFieldManager(G4FieldManager * fM) {
00144 
00145   if (fM!=0) {
00146     fM->SetDetectorField(theField.get());
00147     FieldStepper * theStepper = new FieldStepper(theField->fieldEquation());
00148     theStepper->select(stepper);
00149     G4ChordFinder * CF = new G4ChordFinder(theField.get(),minStep,theStepper);
00150     CF->SetDeltaChord(dChord);
00151     fM->SetChordFinder(CF);
00152     fM->SetDeltaOneStep(dOneStep);
00153     fM->SetDeltaIntersection(dIntersection);
00154     if (dIntersectionAndOneStep != -1.) 
00155       fM->SetAccuraciesWithDeltaOneStep(dIntersectionAndOneStep);
00156   }
00157   if (fChordFinderMonopole == 0) {
00158     G4MonopoleEquation* fMonopoleEquation = new G4MonopoleEquation(theField.get());
00159     G4MagIntegratorStepper* theStepper = new G4ClassicalRK4(fMonopoleEquation,8);
00160     fChordFinderMonopole = new G4ChordFinder(theField.get(),minStep,theStepper);
00161     fChordFinderMonopole->SetDeltaChord(dChord);
00162   }
00163 }
00164 
00165 void FieldBuilder::configurePropagatorInField(G4PropagatorInField * fP) {
00166   if (fP==0) return;
00167   fP->SetMaxLoopCount(int(maxLoopCount));
00168   fP->SetMinimumEpsilonStep(minEpsilonStep);
00169   fP->SetMaximumEpsilonStep(maxEpsilonStep);
00170   fP->SetVerboseLevel(0);
00171   return;
00172 }
00173