CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_8_patch3/src/SimG4CMS/Calo/src/CaloTrkProcessing.cc

Go to the documentation of this file.
00001 #include "SimG4Core/Notification/interface/BeginOfEvent.h"
00002 #include "SimG4Core/Notification/interface/TrackWithHistory.h"
00003 #include "SimG4Core/Notification/interface/TrackInformation.h"
00004 
00005 #include "SimG4CMS/Calo/interface/CaloTrkProcessing.h"
00006 
00007 #include "SimG4Core/Application/interface/SimTrackManager.h"
00008 #include "DetectorDescription/Core/interface/DDFilter.h"
00009 #include "DetectorDescription/Core/interface/DDFilteredView.h"
00010 #include "DetectorDescription/Core/interface/DDSolid.h"
00011 #include "DetectorDescription/Core/interface/DDValue.h"
00012 #include "FWCore/Utilities/interface/Exception.h"
00013 
00014 #include "G4EventManager.hh"
00015 #include "G4LogicalVolumeStore.hh"
00016 #include "G4LogicalVolume.hh"
00017 #include "G4Step.hh"
00018 #include "G4Track.hh"
00019 
00020 //#define DebugLog
00021 
00022 CaloTrkProcessing::CaloTrkProcessing(G4String name, 
00023                                      const DDCompactView & cpv,
00024                                      SensitiveDetectorCatalog & clg, 
00025                                      edm::ParameterSet const & p,
00026                                      const SimTrackManager* manager) : 
00027   SensitiveCaloDetector(name, cpv, clg, p), lastTrackID(-1),
00028   m_trackManager(manager) {  
00029 
00030   //Initialise the parameter set
00031   edm::ParameterSet m_p = p.getParameter<edm::ParameterSet>("CaloTrkProcessing");
00032   testBeam   = m_p.getParameter<bool>("TestBeam");
00033   eMin       = m_p.getParameter<double>("EminTrack")*MeV;
00034   putHistory = m_p.getParameter<bool>("PutHistory");
00035 
00036   edm::LogInfo("CaloSim") << "CaloTrkProcessing: Initailised with TestBeam = " 
00037                           << testBeam << " Emin = " << eMin << " MeV and"
00038                           << " History flag " << putHistory;
00039 
00040   //Get the names 
00041   G4String attribute = "ReadOutName"; 
00042   DDSpecificsFilter filter;
00043   DDValue           ddv(attribute,name,0);
00044   filter.setCriteria(ddv,DDSpecificsFilter::equals);
00045   DDFilteredView fv(cpv);
00046   fv.addFilter(filter);
00047   fv.firstChild();
00048   DDsvalues_type sv(fv.mergedSpecifics());
00049 
00050   G4String value                     = "Calorimeter";
00051   std::vector<std::string> caloNames = getNames (value, sv);
00052 #ifdef DebugLog
00053   LogDebug("CaloSim") << "CaloTrkProcessing: Names for " << value << ":";
00054   for (unsigned int i=0; i<caloNames.size(); i++)
00055     LogDebug("CaloSim") << " (" << i << ") " << caloNames[i];
00056 #endif
00057 
00058   value                                = "Levels";
00059   std::vector<double>      levels      = getNumbers (value, sv);
00060 #ifdef DebugLog
00061   LogDebug("CaloSim") << "CaloTrkProcessing: Names for " << value << ":";
00062   for (unsigned int i=0; i<levels.size(); i++)
00063     LogDebug("CaloSim") << " (" << i << ") " << levels[i];
00064 #endif
00065 
00066   value                                = "Neighbours";
00067   std::vector<double>      neighbours  = getNumbers (value, sv);
00068 #ifdef DebugLog
00069   LogDebug("CaloSim") << "CaloTrkProcessing: Names for " << value << ":";
00070   for (unsigned int i=0; i<neighbours.size(); i++)
00071     LogDebug("CaloSim") << " (" << i << ") " << neighbours[i];
00072 #endif
00073 
00074   value                                = "Inside";
00075   std::vector<std::string> insideNames = getNames (value, sv);
00076 #ifdef DebugLog
00077   LogDebug("CaloSim") << "CaloTrkProcessing: Names for " << value << ":";
00078   for (unsigned int i=0; i<insideNames.size(); i++)
00079     LogDebug("CaloSim") << " (" << i << ") " << insideNames[i];
00080 #endif
00081 
00082   value                                = "InsideLevel";
00083   std::vector<double>      insideLevel = getNumbers (value, sv);
00084 #ifdef DebugLog
00085   LogDebug("CaloSim") << "CaloTrkProcessing: Names for " << value << ":";
00086   for (unsigned int i=0; i<insideLevel.size(); i++)
00087     LogDebug("CaloSim") << " (" << i << ") " << insideLevel[i];
00088 #endif
00089 
00090   if (caloNames.size() < neighbours.size()) {
00091     edm::LogError("CaloSim") << "CaloTrkProcessing: # of Calorimeter bins " 
00092                              << caloNames.size() << " does not match with "
00093                              << neighbours.size() << " ==> illegal ";
00094     throw cms::Exception("Unknown", "CaloTrkProcessing")
00095       << "Calorimeter array size does not match with size of neighbours\n";
00096   }
00097 
00098   const G4LogicalVolumeStore * lvs = G4LogicalVolumeStore::GetInstance();
00099   std::vector<G4LogicalVolume *>::const_iterator lvcite;
00100   int istart = 0;
00101   for (unsigned int i=0; i<caloNames.size(); i++) {
00102     G4LogicalVolume* lv     = 0;
00103     G4String         name   = caloNames[i];
00104     int              number = static_cast<int>(neighbours[i]);
00105     for (lvcite = lvs->begin(); lvcite != lvs->end(); lvcite++) {
00106       if ((*lvcite)->GetName() == name) {
00107         lv = (*lvcite);
00108         break;
00109       }
00110     }
00111     if (lv != 0) {
00112      CaloTrkProcessing::Detector detector;
00113      detector.name  = name;
00114      detector.lv    = lv;
00115      detector.level = static_cast<int>(levels[i]);
00116      if (istart+number > (int)(insideNames.size())) {
00117        edm::LogError("CaloSim") << "CaloTrkProcessing: # of InsideNames bins " 
00118                                 << insideNames.size() <<" too few compaerd to "
00119                                 << istart+number << " requested ==> illegal ";
00120        throw cms::Exception("Unknown", "CaloTrkProcessing")
00121          << "InsideNames array size does not match with list of neighbours\n";
00122      }
00123      std::vector<std::string>      inside;
00124      std::vector<G4LogicalVolume*> insideLV;
00125      std::vector<int>              insideLevels;
00126      for (int k = 0; k < number; k++) {
00127        lv   = 0;
00128        name = insideNames[istart+k];
00129        for (lvcite = lvs->begin(); lvcite != lvs->end(); lvcite++) 
00130          if ((*lvcite)->GetName() == name) {
00131            lv = (*lvcite);
00132            break;
00133          }
00134        inside.push_back(name);
00135        insideLV.push_back(lv);
00136        insideLevels.push_back(static_cast<int>(insideLevel[istart+k]));
00137      }
00138      detector.fromDets   = inside;
00139      detector.fromDetL   = insideLV;
00140      detector.fromLevels = insideLevels;
00141      detectors.push_back(detector);
00142     }
00143     istart += number;
00144   }
00145 
00146   edm::LogInfo("CaloSim") << "CaloTrkProcessing: with " << detectors.size()
00147                           << " calorimetric volumes";
00148   for (unsigned int i=0; i<detectors.size(); i++) {
00149     edm::LogInfo("CaloSim") << "CaloTrkProcessing: Calorimeter volume " << i
00150                             << " " << detectors[i].name << " LV "
00151                             << detectors[i].lv << " at level "
00152                             << detectors[i].level << " with "
00153                             << detectors[i].fromDets.size() << " neighbours";
00154     for (unsigned int k=0; k<detectors[i].fromDets.size(); k++) 
00155       edm::LogInfo("CaloSim") << "                   Element " << k << " "
00156                               << detectors[i].fromDets[k] << " LV "
00157                               << detectors[i].fromDetL[k] << " at level "
00158                               << detectors[i].fromLevels[k];
00159   }
00160 }
00161 
00162 CaloTrkProcessing::~CaloTrkProcessing() {
00163   edm::LogInfo("CaloSim") << "CaloTrkProcessing: Deleted";
00164 }
00165 
00166 void CaloTrkProcessing::update(const BeginOfEvent * evt) {
00167   lastTrackID = -1;
00168 }
00169 
00170 void CaloTrkProcessing::update(const G4Step * aStep) {
00171   
00172   // define if you are at the surface of CALO  
00173   
00174   G4Track* theTrack = aStep->GetTrack();   
00175   int      id       = theTrack->GetTrackID();
00176 
00177   TrackInformation* trkInfo = dynamic_cast<TrackInformation*>
00178     (theTrack->GetUserInformation());
00179   
00180   if (trkInfo == 0) {
00181     edm::LogError("CaloSim") << "CaloTrkProcessing: No trk info !!!! abort ";
00182     throw cms::Exception("Unknown", "CaloTrkProcessing")
00183       << "cannot get trkInfo for Track " << id << "\n";
00184   }
00185   
00186   if (testBeam) {
00187     if (trkInfo->getIDonCaloSurface() == 0) {
00188 #ifdef DebugLog
00189       LogDebug("CaloSim") << "CaloTrkProcessing set IDonCaloSurface to " << id 
00190                           << " at step Number "
00191                           << theTrack->GetCurrentStepNumber();
00192 #endif
00193       trkInfo->setIDonCaloSurface(id,0,0,
00194                                   theTrack->GetDefinition()->GetPDGEncoding(),
00195                                   theTrack->GetMomentum().mag());
00196       lastTrackID = id;
00197       if (theTrack->GetKineticEnergy()/MeV > eMin)
00198         trkInfo->putInHistory();
00199     } 
00200   } else {
00201     if (putHistory) {
00202       trkInfo->putInHistory();
00203       //      trkInfo->setAncestor();
00204     }
00205 #ifdef DebugLog
00206     LogDebug("CaloSim") << "CaloTrkProcessing Entered for " << id 
00207                         << " at stepNumber "<< theTrack->GetCurrentStepNumber()
00208                         << " IDonCaloSur.. " << trkInfo->getIDonCaloSurface()
00209                         << " CaloCheck " << trkInfo->caloIDChecked();
00210 #endif
00211     if (trkInfo->getIDonCaloSurface() != 0) {
00212       if (trkInfo->caloIDChecked() == false) {
00213         G4StepPoint*        postStepPoint = aStep->GetPostStepPoint();   
00214         const G4VTouchable* post_touch    = postStepPoint->GetTouchable();
00215 
00216         if (isItInside(post_touch, trkInfo->getIDCaloVolume(),
00217                        trkInfo->getIDLastVolume()) > 0) {
00218           trkInfo->setIDonCaloSurface(0, -1, -1, 0, 0);
00219         } else {
00220           trkInfo->setCaloIDChecked(true);
00221         }
00222       }
00223     } else {
00224       G4StepPoint*        postStepPoint = aStep->GetPostStepPoint();   
00225       const G4VTouchable* post_touch    = postStepPoint->GetTouchable();
00226       int                 ical          = isItCalo(post_touch);
00227       if (ical >= 0) {
00228         G4StepPoint*        preStepPoint = aStep->GetPreStepPoint(); 
00229         const G4VTouchable* pre_touch    = preStepPoint->GetTouchable();
00230         int                 inside       = isItInside(pre_touch, ical, -1);
00231         if (inside >= 0 ||  (theTrack->GetCurrentStepNumber()==1)) {
00232           trkInfo->setIDonCaloSurface(id,ical,inside,
00233                                       theTrack->GetDefinition()->GetPDGEncoding(),
00234                                       theTrack->GetMomentum().mag());
00235           trkInfo->setCaloIDChecked(true);
00236           lastTrackID = id;
00237           if (theTrack->GetKineticEnergy()/MeV > eMin) trkInfo->putInHistory();
00238 #ifdef DebugLog
00239           LogDebug("CaloSim") <<"CaloTrkProcessing: set ID on Calo " << ical
00240                               << " surface (Inside " << inside << ") to " 
00241                               << id << " of a Track with Kinetic Energy " 
00242                               << theTrack->GetKineticEnergy()/MeV << " MeV";
00243 #endif
00244         }
00245       }
00246     }
00247   }
00248 }
00249 
00250 std::vector<std::string> CaloTrkProcessing::getNames(const G4String str,
00251                                                      const DDsvalues_type &sv){
00252 
00253 #ifdef DebugLog
00254   LogDebug("CaloSim") << "CaloTrkProcessing::getNames called for " << str;
00255 #endif
00256   DDValue value(str);
00257   if (DDfetch(&sv,value)) {
00258 #ifdef DebugLog
00259     LogDebug("CaloSim") << value;
00260 #endif
00261     const std::vector<std::string> & fvec = value.strings();
00262     int nval = fvec.size();
00263     if (nval < 1) {
00264 
00265         edm::LogError("CaloSim") << "CaloTrkProcessing: # of " << str 
00266                                  << " bins " << nval << " < 1 ==> illegal ";
00267         throw cms::Exception("Unknown", "CaloTrkProcessing")
00268           << "nval < 2 for array " << str << "\n";
00269     }
00270     
00271     return fvec;
00272   } else {
00273     edm::LogError("CaloSim") << "CaloTrkProcessing: cannot get array " << str ;
00274     throw cms::Exception("Unknown", "CaloTrkProcessing")
00275       << "cannot get array " << str << "\n";
00276   }
00277 }
00278 
00279 std::vector<double> CaloTrkProcessing::getNumbers(const G4String str,
00280                                                   const DDsvalues_type &sv) {
00281 
00282 #ifdef DebugLog
00283   LogDebug("CaloSim") << "CaloTrkProcessing::getNumbers called for " << str;
00284 #endif
00285   DDValue value(str);
00286   if (DDfetch(&sv,value)) {
00287 #ifdef DebugLog
00288     LogDebug("CaloSim") << value;
00289 #endif
00290     const std::vector<double> & fvec = value.doubles();
00291     int nval = fvec.size();
00292     if (nval < 1) {
00293         edm::LogError("CaloSim") << "CaloTrkProcessing: # of " << str 
00294                                  << " bins " << nval << " < 1 ==> illegal ";
00295         throw cms::Exception("Unknown", "CaloTrkProcessing")
00296           << "nval < 2 for array " << str << "\n";
00297     }
00298     
00299     return fvec;
00300   } else {
00301     edm::LogError("CaloSim") << "CaloTrkProcessing: cannot get array " << str ;
00302     throw cms::Exception("Unknown", "CaloTrkProcessing")
00303       << "cannot get array " << str << "\n";
00304   }
00305 }
00306 
00307 int CaloTrkProcessing::isItCalo(const G4VTouchable* touch) {
00308 
00309   int lastLevel = -1;
00310   G4LogicalVolume* lv=0;
00311   for (unsigned int it=0; it < detectors.size(); it++) {
00312     if (lastLevel != detectors[it].level) {
00313       lastLevel = detectors[it].level;
00314       lv        = detLV(touch, lastLevel);
00315 #ifdef DebugLog
00316       std::string  name1 = "Unknown";
00317       if (lv != 0) name1 = lv->GetName(); 
00318       LogDebug("CaloSim") << "CaloTrkProcessing: volume " << name1
00319                           << " at Level " << lastLevel;
00320       int levels = detLevels(touch);
00321       if (levels > 0) {
00322         G4String name2[20]; int copyno2[20];
00323         detectorLevel(touch, levels, copyno2, name2);
00324         for (int i2=0; i2<levels; i2++) 
00325           LogDebug("CaloSim") << " " << i2 << " " << name2[i2] << " " 
00326                               << copyno2[i2];
00327       }
00328 #endif
00329     }
00330     bool ok = (lv   == detectors[it].lv);
00331     if (ok) return it;
00332   }
00333   return -1;
00334 }
00335 
00336 int CaloTrkProcessing::isItInside(const G4VTouchable* touch, int idcal,
00337                                   int idin) {
00338   int lastLevel = -1;
00339   G4LogicalVolume* lv=0;
00340   int id1, id2;
00341   if (idcal < 0) {id1 = 0; id2 = static_cast<int>(detectors.size());}
00342   else           {id1 = idcal; id2 = id1+1;}
00343   for (int it1 = id1; it1 < id2; it1++) {
00344     if (idin < 0) {
00345       for (unsigned int it2 = 0; it2 < detectors[it1].fromDets.size(); it2++) {
00346         if (lastLevel != detectors[it1].fromLevels[it2]) {
00347           lastLevel = detectors[it1].fromLevels[it2];
00348           lv        = detLV(touch,lastLevel);
00349 #ifdef DebugLog
00350           std::string  name1 = "Unknown";
00351           if (lv != 0) name1 = lv->GetName(); 
00352           LogDebug("CaloSim") << "CaloTrkProcessing: volume " << name1
00353                               << " at Level " << lastLevel;
00354           int levels = detLevels(touch);
00355           if (levels > 0) {
00356             G4String name2[20]; int copyno2[20];
00357             detectorLevel(touch, levels, copyno2, name2);
00358             for (int i2=0; i2<levels; i2++) 
00359               LogDebug("CaloSim") << " " << i2 << " " << name2[i2] << " " 
00360                                   << copyno2[i2];
00361           }
00362 #endif
00363         }
00364         bool ok = (lv   == detectors[it1].fromDetL[it2]);
00365         if (ok) return it2;
00366       }
00367     } else {
00368       lastLevel = detectors[it1].fromLevels[idin];
00369       lv        = detLV(touch,lastLevel);
00370 #ifdef DebugLog
00371       std::string  name1 = "Unknown";
00372       if (lv != 0) name1 = lv->GetName(); 
00373       LogDebug("CaloSim") << "CaloTrkProcessing: volume " << name1
00374                           << " at Level " << lastLevel;
00375       int levels = detLevels(touch);
00376       if (levels > 0) {
00377         G4String name2[20]; int copyno2[20];
00378         detectorLevel(touch, levels, copyno2, name2);
00379         for (int i2=0; i2<levels; i2++) 
00380           LogDebug("CaloSim") << " " << i2 << " " << name2[i2] << " " 
00381                               << copyno2[i2];
00382       }
00383 #endif
00384       bool ok = (lv   == detectors[it1].fromDetL[idin]);
00385       if (ok) return idin;
00386     }
00387   }
00388   return -1;
00389 }
00390 
00391 int CaloTrkProcessing::detLevels(const G4VTouchable* touch) const {
00392 
00393   //Return number of levels
00394   if (touch) 
00395     return ((touch->GetHistoryDepth())+1);
00396   else
00397     return 0;
00398 }
00399 
00400 G4LogicalVolume* CaloTrkProcessing::detLV(const G4VTouchable* touch,
00401                                           int currentlevel) const {
00402 
00403   G4LogicalVolume* lv=0;
00404   if (touch) {
00405     int level = ((touch->GetHistoryDepth())+1);
00406     if (level > 0 && level >= currentlevel) {
00407       int ii = level - currentlevel; 
00408       lv     = touch->GetVolume(ii)->GetLogicalVolume();
00409       return lv;
00410     } 
00411   }
00412   return lv;
00413 }
00414 
00415 void CaloTrkProcessing::detectorLevel(const G4VTouchable* touch, int& level,
00416                                       int* copyno, G4String* name) const {
00417 
00418   static const std::string unknown("Unknown");
00419   //Get name and copy numbers
00420   if (level > 0) {
00421     for (int ii = 0; ii < level; ii++) {
00422       int i      = level - ii - 1;
00423       G4VPhysicalVolume* pv = touch->GetVolume(i);
00424       if (pv != 0) 
00425         name[ii] = pv->GetName();
00426       else
00427         name[ii] = unknown;
00428       copyno[ii] = touch->GetReplicaNumber(i);
00429     }
00430   }
00431 }