CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/SimG4CMS/Tracker/src/TouchableToHistory.cc

Go to the documentation of this file.
00001 #include "SimG4CMS/Tracker/interface/TouchableToHistory.h"
00002 #include "DetectorDescription/Core/interface/DDExpandedView.h"
00003 #include "DetectorDescription/Core/interface/DDFilteredView.h"
00004 #include "DetectorDescription/Core/interface/DDCompactView.h"
00005 #include "Geometry/TrackerNumberingBuilder/interface/GeometricDet.h"
00006 
00007 
00008 #include "G4TransportationManager.hh" 
00009 #include "G4Navigator.hh" 
00010 #include "G4VTouchable.hh"
00011 #include "G4TouchableHistory.hh"
00012 
00013 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00014 
00015 //#define DEBUG
00016 
00017 void TouchableToHistory::buildAll(){
00018   if (alreadySet) return;
00019   alreadySet = true;
00020 
00021   G4Navigator* theStdNavigator = G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking();
00022   G4Navigator theNavigator;
00023   theNavigator.SetWorldVolume(theStdNavigator->GetWorldVolume());
00024 
00025 
00026   std::vector<const GeometricDet*> allSensitiveDets;
00027   myGeomDet->deepComponents(allSensitiveDets);
00028   edm::LogInfo("TrackerSimInfoNumbering")<<" TouchableTo History: got "<<allSensitiveDets.size()<<" sensitive detectors from TrackerMapDDDtoID.";
00029 
00030   for ( std::vector<const GeometricDet*>::const_iterator it = allSensitiveDets.begin(); 
00031         it != allSensitiveDets.end(); 
00032         ++it)
00033     {
00034       DDTranslation const & t = (*it)->translation(); 
00035       theNavigator.LocateGlobalPointAndSetup(G4ThreeVector(t.x(),t.y(),t.z()));
00036       G4TouchableHistory * hist = theNavigator.CreateTouchableHistory(); 
00037       TouchableToHistory::Nav_Story st =  touchableToNavStory(hist);
00038 
00039 #ifdef DEBUG    
00040     u_int32_t oldsize = myDirectMap.size();
00041 #endif
00042 
00043     myMap[st] = nav_type((*it)->navType().begin(),(*it)->navType().end());
00044     myDirectMap[st] = (*it)->geographicalID();
00045 
00046     /*
00047 #ifdef DEBUG    
00048     LogDebug("TrackerSimDebugNumbering")<< " INSERTING "<<view.logicalPart().name()<<" "<<t<<" "<<hist->GetVolume()->GetLogicalVolume()->GetName();
00049     LogDebug("TrackerSimDebugNumbering")<<" Sensitive: "<<hist->GetVolume()->GetLogicalVolume()->GetSensitiveDetector()<<std::endl;
00050     LogDebug("TrackerSimDebugNumbering")<<"Now size is "<<myDirectMap.size()<<std::endl;
00051     if (oldsize == myDirectMap.size())
00052       edm::LogError("TrackerSimInfoNumbering")<< "Touchable to History Error!!!!";
00053     dumpG4VPV(hist);
00054 #endif
00055     */
00056     delete hist;
00057 
00058   }
00059   edm::LogInfo("TrackerSimInfoNumbering")<<" TouchableToHistory: mapped "<<myDirectMap.size()<<" detectors to G4.";
00060 
00061   if (myDirectMap.size() != allSensitiveDets.size()){
00062     edm::LogError("TrackerSimInfoNumbering")<<" ERROR: DDD sensitive detectors do not match Geant4 ones.";
00063     //FIXME use throw
00064     abort();
00065   }
00066 
00067 
00068 }
00069 
00070 DDFilteredView& TouchableToHistory::getFilteredView(const G4VTouchable& t, DDFilteredView& f){
00071   if (alreadySet == false) 
00072     buildAll();
00073   f.goTo(myMap[touchableToNavStory(&t)]);
00074   return f;
00075 }
00076 TouchableToHistory::nav_type TouchableToHistory::getNavType(const G4VTouchable& t){
00077   if (alreadySet == false) 
00078     edm::LogError("TrackerSimInfoNumbering")<<" NOT READY ";
00079   return myMap[touchableToNavStory(&t)];
00080 }
00081 
00082 TouchableToHistory::Nav_Story TouchableToHistory::getNavStory(DDFilteredView& i){
00083   if (alreadySet == false) buildAll();
00084   DDTranslation t = i.translation();
00085 
00086   G4Navigator* theStdNavigator = G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking();
00087   G4Navigator  theNavigator;
00088   theNavigator.SetWorldVolume(theStdNavigator->GetWorldVolume());
00089 
00090   theNavigator.LocateGlobalPointAndSetup(G4ThreeVector(t.x(),t.y(),t.z()));
00091   G4TouchableHistory* hist = theNavigator.CreateTouchableHistory(); 
00092   TouchableToHistory::Nav_Story temp = touchableToNavStory(hist);
00093   delete hist;
00094   return (temp);
00095 }
00096 
00097 TouchableToHistory::Nav_Story TouchableToHistory::touchableToNavStory(const G4VTouchable *v) {
00098   static G4String tobinactive("TOBInactive");
00099   Nav_Story temp;
00100 #ifdef DEBUG    
00101   std::vector<int> debugint;
00102   std::vector<std::string> debugstring;
00103 #endif
00104   int levels = v->GetHistoryDepth();
00105   
00106   for (int k=0; k<=levels; k++){
00107     if (v->GetVolume(k)->GetLogicalVolume()->GetName() != tobinactive) {
00108       temp.push_back(
00109                      std::pair<int,std::string>
00110                      (v->GetVolume(k)->GetCopyNo(),
00111                       v->GetVolume(k)->GetLogicalVolume()->GetName()));
00112 #ifdef DEBUG    
00113       debugint.push_back(v->GetVolume(k)->GetCopyNo());
00114       debugstring.push_back(v->GetVolume(k)->GetLogicalVolume()->GetName());
00115 #endif
00116     }
00117   }
00118 #ifdef DEBUG    
00119   // LogDebug("TrackerSimDebugNumbering")<<" G4 TouchableToHistory "<< debugint;
00120   for(u_int32_t jj=0;jj<debugstring.size();jj++)LogDebug("TrackerSimDebugNumbering")<<" "<<debugstring[jj];
00121 #endif
00122   return temp;
00123 }
00124 
00125 TouchableToHistory::nav_type TouchableToHistory::touchableToNavType(const G4VTouchable* v) {
00126   if (alreadySet == false) 
00127     buildAll();
00128 
00129   dumpG4VPV(v);
00130 
00131   return   myMap[touchableToNavStory(v)];
00132 }
00133 
00134 int TouchableToHistory::touchableToInt(const G4VTouchable* v) {
00135   if (alreadySet == false) 
00136     buildAll();
00137 
00138   dumpG4VPV(v);
00139 
00140   LogDebug("TrackerSimDebugNumbering")<<" Returning: "<< myDirectMap[touchableToNavStory(v)]<<std::endl;
00141 
00142   return   myDirectMap[touchableToNavStory(v)];
00143 }
00144 
00145 void TouchableToHistory::dumpG4VPV(const G4VTouchable* v){
00146   int levels = v->GetHistoryDepth();
00147   
00148   LogDebug("TrackerSimDebugNumbering")<<" NAME : "<<v->GetVolume()->GetLogicalVolume()->GetName();
00149   for (int k=0; k<=levels; k++){
00150     LogDebug("TrackerSimInfoNumbering") <<" Hist: "<< v->GetVolume(k)->GetLogicalVolume()->GetName()<<
00151       " Copy "<< v->GetVolume(k)->GetCopyNo();
00152   }
00153 }
00154