CMS 3D CMS Logo

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/Notification/interface/SimG4Exception.h"
00013 
00014 /*
00015 #include "DetectorDescription/Parser/interface/DDLConfiguration.h"
00016 #include "DetectorDescription/Base/interface/DDException.h"
00017 #include "DetectorDescription/Algorithm/src/AlgoInit.h"
00018 */
00019 
00020 #include "G4Mag_UsualEqRhs.hh"
00021 #include "G4PropagatorInField.hh"
00022 #include "G4FieldManager.hh"
00023 #include "G4TransportationManager.hh"
00024 #include "G4ChordFinder.hh"
00025 #include "G4UniformMagField.hh"
00026 
00027 #include "SimG4Core/MagneticField/interface/LocalFieldManager.h"
00028 
00029 #include "G4LogicalVolumeStore.hh"
00030 
00031 using namespace sim;
00032 
00033 FieldBuilder::FieldBuilder(const MagneticField * f, 
00034                            const edm::ParameterSet & p) 
00035    : theField( new Field(f,p.getParameter<double>("delta"))), 
00036      theFieldEquation(new G4Mag_UsualEqRhs(theField.get())),
00037      theTopVolume(0),
00038      fieldValue(0.), minStep(0.), dChord(0.), dOneStep(0.),
00039      dIntersection(0.), dIntersectionAndOneStep(0.), 
00040      maxLoopCount(0), minEpsilonStep(0.), maxEpsilonStep(0.), 
00041      thePSet(p)
00042 {
00043     theField->fieldEquation(theFieldEquation);
00044 }
00045 
00046 /*
00047 void FieldBuilder::readFieldParameters(DDLogicalPart lp,
00048                                        const std::string& keywordField)
00049 {
00050     int tmp;
00051     tmp = map_.toString(keywordField,lp,fieldType);
00052     tmp = map_.toDouble("FieldValue",lp,fieldValue);
00053     tmp = map_.toString("Stepper",lp,stepper);
00054     tmp = map_.toDouble("MinStep",lp,minStep);
00055     tmp = map_.toDouble("DeltaChord",lp,dChord);
00056     tmp = map_.toDouble("DeltaOneStep",lp,dOneStep);
00057     tmp = map_.toDouble("DeltaIntersection",lp,dIntersection);
00058     tmp = map_.toDouble("DeltaIntersectionAndOneStep",lp,dIntersectionAndOneStep);
00059     tmp = map_.toDouble("MaximumLoopCount",lp,maxLoopCount);
00060     tmp = map_.toDouble("MinimumEpsilonStep",lp,minEpsilonStep);
00061     tmp = map_.toDouble("MaximumEpsilonStep",lp,maxEpsilonStep);
00062     LogDebug("MagneticField") << "FieldBuilder " << lp.name() << " " 
00063                               << fieldType << " " << fieldValue << " " 
00064                               << stepper << " " << minStep << " " << dChord 
00065                               << " " << dOneStep << " " << dIntersection << " "
00066                               << dIntersectionAndOneStep << " "
00067                               << maxLoopCount << " " << minEpsilonStep << " " 
00068                               << maxEpsilonStep;
00069     return;
00070 }
00071 */
00072 
00073 void FieldBuilder::build( G4FieldManager* fM, G4PropagatorInField* fP)
00074 {
00075     
00076     edm::ParameterSet thePSetForGMFM =
00077        thePSet.getParameter<edm::ParameterSet>("ConfGlobalMFM");
00078 
00079     std::string volName = thePSetForGMFM.getParameter< std::string >("Volume");
00080 
00081     edm::ParameterSet volPSet =
00082        thePSetForGMFM.getParameter< edm::ParameterSet >( volName );
00083     
00084     configureForVolume( volName, volPSet, fM, fP );
00085     
00086     // configure( "MagneticFieldType", fM, fP ) ;
00087 
00088     if ( thePSet.getParameter<bool>("UseLocalMagFieldManager") )
00089     {
00090 
00091        edm::ParameterSet defpset ;
00092        edm::ParameterSet thePSetForLMFM = 
00093           thePSet.getUntrackedParameter<edm::ParameterSet>("ConfLocalMFM", defpset);
00094        //
00095        // Patology !!! LocalFM requested but configuration not given ! 
00096        // In principal, need to throw an exception
00097        //
00098        if ( thePSetForLMFM == defpset )
00099        {
00100           std::cout << " Patology ! Local Mag.Field Manager requested but config not given !\n";
00101           return ;
00102        }
00103        
00104        std::vector<std::string> ListOfVolumes = 
00105           thePSetForLMFM.getParameter< std::vector<std::string> >("ListOfVolumes");
00106           
00107        // creating Local Mag.Field Manager
00108        for (unsigned int i = 0; i < ListOfVolumes.size(); ++ i )
00109        {
00110           volPSet = thePSetForLMFM.getParameter< edm::ParameterSet >(ListOfVolumes[i]);
00111           G4FieldManager* fAltM = new G4FieldManager() ;
00112           configureForVolume( ListOfVolumes[i], volPSet, fAltM ) ;
00113           //configureLocalFM( ListOfVolumes[i], fAltM ) ;
00114           LocalFieldManager* fLM = new LocalFieldManager( theField.get(), fM, fAltM ) ;
00115           fLM->SetVerbosity(thePSet.getUntrackedParameter<bool>("Verbosity",false));
00116           theTopVolume->SetFieldManager( fLM, true ) ;
00117        }
00118     }
00119     return ;
00120  
00121 }
00122 
00123 /*
00124 void FieldBuilder::configure(const std::string& keywordField,
00125                              G4FieldManager * fM,
00126                              G4PropagatorInField * fP)
00127 {
00128   G4LogicalVolumeToDDLogicalPartMap::Vector vec = map_.all(keywordField);
00129   for (G4LogicalVolumeToDDLogicalPartMap::Vector::iterator tit = vec.begin();
00130        tit != vec.end(); tit++) {
00131     theTopVolume = (*tit).first;
00132     readFieldParameters((*tit).second,keywordField);    
00133     if (fM!=0) configureFieldManager(fM);
00134     if (fP!=0) configurePropagatorInField(fP);  
00135   }
00136   return ;
00137 }
00138 */
00139 
00140 void FieldBuilder::configureForVolume( const std::string& volName,
00141                                        edm::ParameterSet& volPSet,
00142                                        G4FieldManager * fM,
00143                                        G4PropagatorInField * fP )
00144 {
00145 
00146    G4LogicalVolumeStore* theStore = G4LogicalVolumeStore::GetInstance();
00147    for (unsigned int i=0; i<(*theStore).size(); ++i )
00148    {
00149       std::string curVolName = ((*theStore)[i])->GetName();
00150       if ( curVolName == volName )
00151       {
00152          theTopVolume = (*theStore)[i] ;
00153       }
00154    }
00155 
00156    fieldType     = volPSet.getParameter<std::string>("Type") ;
00157    stepper       = volPSet.getParameter<std::string>("Stepper") ;
00158    edm::ParameterSet stpPSet = 
00159       volPSet.getParameter<edm::ParameterSet>(stepper) ;
00160    minStep       = stpPSet.getParameter<double>("MinStep") ;
00161    dChord        = stpPSet.getParameter<double>("DeltaChord") ;
00162    dOneStep      = stpPSet.getParameter<double>("DeltaOneStep") ;
00163    dIntersection = stpPSet.getParameter<double>("DeltaIntersection") ;
00164    dIntersectionAndOneStep = stpPSet.getUntrackedParameter<double>("DeltaIntersectionAndOneStep",-1.);
00165    maxLoopCount = stpPSet.getUntrackedParameter<double>("MaximumLoopCounts",1000);
00166    minEpsilonStep = stpPSet.getUntrackedParameter<double>("MinimumEpsilonStep",0.00001);
00167    maxEpsilonStep = stpPSet.getUntrackedParameter<double>("MaximumEpsilonStep",0.01);
00168    
00169    if (fM!=0) configureFieldManager(fM);
00170    if (fP!=0) configurePropagatorInField(fP);   
00171 
00172    return;
00173 
00174 }
00175 
00176 void FieldBuilder::configureFieldManager(G4FieldManager * fM)
00177 {
00178     if (fM==0) return;
00179     fM->SetDetectorField(theField.get());
00180     FieldStepper * theStepper = new FieldStepper(theField->fieldEquation());
00181     theStepper->select(stepper);
00182     G4ChordFinder * CF = new G4ChordFinder(theField.get(),minStep,theStepper);
00183     CF->SetDeltaChord(dChord);
00184     fM->SetChordFinder(CF);
00185     fM->SetDeltaOneStep(dOneStep);
00186     fM->SetDeltaIntersection(dIntersection);
00187     if (dIntersectionAndOneStep != -1.) 
00188         fM->SetAccuraciesWithDeltaOneStep(dIntersectionAndOneStep);
00189     return;
00190 }
00191 
00192 void FieldBuilder::configurePropagatorInField(G4PropagatorInField * fP)
00193 {
00194     if (fP==0) return;
00195     fP->SetMaxLoopCount(int(maxLoopCount));
00196     fP->SetMinimumEpsilonStep(minEpsilonStep);
00197     fP->SetMaximumEpsilonStep(maxEpsilonStep);
00198     fP->SetVerboseLevel(0);
00199     return;
00200 }
00201 
00202 G4LogicalVolume * FieldBuilder::fieldTopVolume() { return theTopVolume; }

Generated on Tue Jun 9 17:47:05 2009 for CMSSW by  doxygen 1.5.4