CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10_patch2/src/SimG4CMS/Muon/src/MuonSensitiveDetector.cc

Go to the documentation of this file.
00001 #include "SimG4CMS/Muon/interface/MuonSensitiveDetector.h"
00002 #include "SimG4CMS/Muon/interface/MuonSlaveSD.h"
00003 #include "SimG4CMS/Muon//interface/MuonEndcapFrameRotation.h"
00004 #include "SimG4CMS/Muon/interface/MuonRpcFrameRotation.h"
00005 #include "Geometry/MuonNumbering/interface/MuonSubDetector.h"
00006 
00007 #include "DataFormats/GeometryVector/interface/LocalVector.h"
00008 
00009 #include "SimG4CMS/Muon/interface/SimHitPrinter.h"
00010 #include "SimDataFormats/TrackingHit/interface/UpdatablePSimHit.h"
00011 
00012 #include "SimG4CMS/Muon/interface/MuonG4Numbering.h"
00013 #include "Geometry/MuonNumbering/interface/MuonBaseNumber.h"
00014 #include "Geometry/MuonNumbering/interface/MuonSimHitNumberingScheme.h"
00015 
00016 #include "SimG4Core/Notification/interface/TrackInformation.h"
00017 #include "SimG4Core/Notification/interface/G4TrackToParticleID.h"
00018 #include "SimG4Core/Physics/interface/G4ProcessTypeEnumerator.h"
00019 
00020 #include "G4SDManager.hh"
00021 #include "G4VProcess.hh"
00022 #include "G4EventManager.hh"
00023 
00024 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00025 
00026 #include <iostream>
00027 
00028 
00029 MuonSensitiveDetector::MuonSensitiveDetector(std::string name, 
00030                                              const DDCompactView & cpv,
00031                                              SensitiveDetectorCatalog & clg,
00032                                              edm::ParameterSet const & p,
00033                                              const SimTrackManager* manager) 
00034   : SensitiveTkDetector(name, cpv, clg, p),
00035     thePV(0), theHit(0), theDetUnitId(0), theTrackID(0), theManager(manager)
00036 {
00037   edm::ParameterSet m_MuonSD = p.getParameter<edm::ParameterSet>("MuonSD");
00038   STenergyPersistentCut = m_MuonSD.getParameter<double>("EnergyThresholdForPersistency");//Default 1. GeV
00039   STallMuonsPersistent = m_MuonSD.getParameter<bool>("AllMuonsPersistent");
00040   printHits  = m_MuonSD.getParameter<bool>("PrintHits");
00041   
00042   //
00043   // Here simply create 1 MuonSlaveSD for the moment
00044   //  
00045   
00046   LogDebug("MuonSimDebug") << "create MuonSubDetector "<<name<<std::endl;
00047 
00048   detector = new MuonSubDetector(name);
00049 
00050   LogDebug("MuonSimDebug") << "create MuonFrameRotation"<<std::endl;
00051 
00052  if (detector->isEndcap()) {
00053    //    cout << "MuonFrameRotation create MuonEndcapFrameRotation"<<endl;
00054     theRotation=new MuonEndcapFrameRotation();
00055   } else if (detector->isRpc()) {
00056     //    cout << "MuonFrameRotation create MuonRpcFrameRotation"<<endl;
00057     theRotation=new MuonRpcFrameRotation( cpv );
00058   }  else {
00059     theRotation = 0;
00060   }
00061   LogDebug("MuonSimDebug") << "create MuonSlaveSD"<<std::endl;
00062   slaveMuon  = new MuonSlaveSD(detector,theManager);
00063   LogDebug("MuonSimDebug") << "create MuonSimHitNumberingScheme"<<std::endl;
00064   numbering  = new MuonSimHitNumberingScheme(detector, cpv);
00065   g4numbering = new MuonG4Numbering(cpv);
00066   
00067 
00068   //
00069   // Now attach the right detectors (LogicalVolumes) to me
00070   //
00071   std::vector<std::string>  lvNames = clg.logicalNames(name);
00072   this->Register();
00073   for (std::vector<std::string>::iterator it = lvNames.begin();  it != lvNames.end(); it++){
00074     LogDebug("MuonSimDebug") << name << " MuonSensitiveDetector:: attaching SD to LV " << *it << std::endl;
00075     this->AssignSD(*it);
00076   }
00077 
00078   if (printHits) {
00079     thePrinter = new SimHitPrinter("HitPositionOSCAR.dat");
00080   }
00081 
00082 
00083     LogDebug("MuonSimDebug") << "  EnergyThresholdForPersistency " << STenergyPersistentCut << " AllMuonsPersistent " <<  STallMuonsPersistent << std::endl;
00084     
00085     theG4ProcessTypeEnumerator = new G4ProcessTypeEnumerator;
00086     myG4TrackToParticleID = new G4TrackToParticleID;
00087 
00088 }
00089 
00090 
00091 MuonSensitiveDetector::~MuonSensitiveDetector() { 
00092   delete g4numbering;
00093   delete numbering;
00094   delete slaveMuon;
00095   delete theRotation;
00096   delete detector;
00097 
00098   delete theG4ProcessTypeEnumerator;
00099   
00100   delete myG4TrackToParticleID;
00101 }
00102 
00103 void MuonSensitiveDetector::update(const BeginOfEvent * i){
00104   clearHits();
00105 
00106   //----- Initialize variables to check if two steps belong to same hit
00107   thePV = 0;
00108   theDetUnitId = 0;
00109   theTrackID = 0;
00110 
00111 }
00112 
00113 void MuonSensitiveDetector::update(const  ::EndOfEvent * ev)
00114 {
00115   //slaveMuon->renumbering(theManager);
00116 }
00117 
00118 
00119 void MuonSensitiveDetector::clearHits()
00120 {
00121   LogDebug("MuonSimDebug") << "MuonSensitiveDetector::clearHits"<<std::endl;
00122   slaveMuon->Initialize();
00123 }
00124 
00125 bool MuonSensitiveDetector::ProcessHits(G4Step * aStep, G4TouchableHistory * ROhist)
00126 {
00127   LogDebug("MuonSimDebug") <<" MuonSensitiveDetector::ProcessHits "<<InitialStepPosition(aStep,WorldCoordinates)<<std::endl;
00128 
00129  // TimeMe t1( theHitTimer, false);
00130 
00131   if (aStep->GetTotalEnergyDeposit()>0.){
00132     // do not count neutrals that are killed by User Limits MinEKine
00133     if( aStep->GetTrack()->GetDynamicParticle()->GetCharge() != 0 ){
00134   
00135       if (newHit(aStep)) {
00136         saveHit();
00137         createHit(aStep);
00138       } else {
00139         updateHit(aStep);
00140       }    
00141       return true;
00142     } else {
00143       storeVolumeAndTrack(aStep);
00144       return false;
00145     }
00146   }
00147   return false;
00148 }
00149 
00150 uint32_t MuonSensitiveDetector::setDetUnitId(G4Step * aStep)
00151 { 
00152   //  G4VPhysicalVolume * pv = aStep->GetPreStepPoint()->GetPhysicalVolume();
00153   MuonBaseNumber num = g4numbering->PhysicalVolumeToBaseNumber(aStep);
00154   return numbering->baseNumberToUnitNumber(num);
00155 }
00156 
00157 
00158 Local3DPoint MuonSensitiveDetector::toOrcaRef(Local3DPoint in ,G4Step * s){
00159   if (theRotation !=0 ) {
00160     return theRotation->transformPoint(in,s);
00161   }
00162   return (in);
00163 }
00164 
00165 Local3DPoint MuonSensitiveDetector::toOrcaUnits(Local3DPoint in){
00166   return Local3DPoint(in.x()/cm,in.y()/cm,in.z()/cm);
00167 }
00168 
00169 Global3DPoint MuonSensitiveDetector::toOrcaUnits(Global3DPoint in){
00170   return Global3DPoint(in.x()/cm,in.y()/cm,in.z()/cm);
00171 }
00172 
00173 void MuonSensitiveDetector::storeVolumeAndTrack(G4Step * aStep)
00174 {
00175   G4VPhysicalVolume* pv = aStep->GetPreStepPoint()->GetPhysicalVolume();
00176   G4Track * t  = aStep->GetTrack();
00177   thePV=pv;
00178   theTrackID=t->GetTrackID();
00179 }
00180 
00181 bool MuonSensitiveDetector::newHit(G4Step * aStep){
00182   
00183   G4VPhysicalVolume* pv = aStep->GetPreStepPoint()->GetPhysicalVolume();
00184   G4Track * t  = aStep->GetTrack();
00185   uint32_t currentUnitId=setDetUnitId(aStep);
00186   unsigned int currentTrackID=t->GetTrackID();
00187   //unsigned int currentEventID=G4EventManager::GetEventManager()->GetConstCurrentEvent()->GetEventID();
00188   bool changed=((pv!=thePV) || 
00189                 (currentUnitId!=theDetUnitId) || 
00190                 (currentTrackID!=theTrackID));
00191   return changed;
00192 }
00193 
00194 void MuonSensitiveDetector::createHit(G4Step * aStep){
00195 
00196   G4Track * theTrack  = aStep->GetTrack(); 
00197 
00198   Local3DPoint theEntryPoint;
00199   Local3DPoint theExitPoint;
00200 
00201   if (detector->isBarrel()) {
00202     theEntryPoint= toOrcaUnits(toOrcaRef(InitialStepPositionVsParent(aStep,1),aStep)); // 1 level up
00203     theExitPoint= toOrcaUnits(toOrcaRef(FinalStepPositionVsParent(aStep,1),aStep));  
00204   } else if (detector->isEndcap()) {
00205     // save local z at current level
00206     theEntryPoint= toOrcaUnits(toOrcaRef(InitialStepPosition(aStep,LocalCoordinates),aStep));
00207     theExitPoint= toOrcaUnits(toOrcaRef(FinalStepPosition(aStep,LocalCoordinates),aStep));
00208     float zentry = theEntryPoint.z();
00209     float zexit = theExitPoint.z();
00210     Local3DPoint tempEntryPoint= toOrcaUnits(toOrcaRef(InitialStepPositionVsParent(aStep,4),aStep)); // 4 levels up
00211     Local3DPoint tempExitPoint= toOrcaUnits(toOrcaRef(FinalStepPositionVsParent(aStep,4),aStep));
00212     // reset local z from z wrt deep-parent volume to z wrt low-level volume
00213     theEntryPoint = Local3DPoint( tempEntryPoint.x(), tempEntryPoint.y(), zentry );
00214     theExitPoint  = Local3DPoint( tempExitPoint.x(), tempExitPoint.y(), zexit );
00215   } else {
00216     theEntryPoint= toOrcaUnits(toOrcaRef(InitialStepPosition(aStep,LocalCoordinates),aStep));
00217     theExitPoint= toOrcaUnits(toOrcaRef(FinalStepPosition(aStep,LocalCoordinates),aStep)); 
00218   }
00219 
00220   float thePabs             = aStep->GetPreStepPoint()->GetMomentum().mag()/GeV;
00221   float theTof              = aStep->GetPreStepPoint()->GetGlobalTime()/nanosecond;
00222   float theEnergyLoss       = aStep->GetTotalEnergyDeposit()/GeV;
00223   //  int theParticleType     = theTrack->GetDefinition()->GetPDGEncoding();
00224   int theParticleType     = myG4TrackToParticleID->particleID(theTrack);
00225   G4ThreeVector gmd  = aStep->GetPreStepPoint()->GetMomentumDirection();
00226   G4ThreeVector lmd = ((G4TouchableHistory *)(aStep->GetPreStepPoint()->GetTouchable()))->GetHistory()
00227     ->GetTopTransform().TransformAxis(gmd);
00228   Local3DPoint lnmd = toOrcaRef(ConvertToLocal3DPoint(lmd),aStep);
00229   float theThetaAtEntry = lnmd.theta();
00230   float thePhiAtEntry = lnmd.phi();
00231 
00232   storeVolumeAndTrack( aStep );
00233   theDetUnitId              = setDetUnitId(aStep);
00234 
00235   Global3DPoint theGlobalPos;
00236   if (printHits) {   
00237     Local3DPoint theGlobalHelp = InitialStepPosition(aStep,WorldCoordinates);
00238     theGlobalEntry = toOrcaUnits(Global3DPoint (theGlobalHelp.x(),theGlobalHelp.y(),theGlobalHelp.z()));
00239 
00240     G4StepPoint * preStepPoint = aStep->GetPreStepPoint();
00241     G4TouchableHistory * theTouchable=(G4TouchableHistory *)
00242                                       (preStepPoint->GetTouchable());
00243     theGlobalHelp=ConvertToLocal3DPoint(theTouchable->GetTranslation());
00244     theGlobalPos = toOrcaUnits(Global3DPoint (theGlobalHelp.x(),theGlobalHelp.y(),theGlobalHelp.z()));
00245     //    const G4RotationMatrix * theGlobalRot = theTouchable->GetRotation();
00246   }
00247   
00248   LogDebug("MuonSimDebug") << "MuonSensitiveDetector::createHit UpdatablePSimHit"<<std::endl;
00249 
00250   theHit = new UpdatablePSimHit(theEntryPoint,theExitPoint,thePabs,theTof,
00251                   theEnergyLoss,theParticleType,theDetUnitId,
00252                   theTrackID,theThetaAtEntry,thePhiAtEntry,
00253                   theG4ProcessTypeEnumerator->processId(theTrack->GetCreatorProcess()));
00254 
00255 
00256   LogDebug("MuonSimDebug") <<"=== NEW ==================> ELOSS   = "<<theEnergyLoss<<" "
00257        <<thePV->GetLogicalVolume()->GetName()<<std::endl;
00258   const G4VProcess* p = aStep->GetPostStepPoint()->GetProcessDefinedStep();
00259   const G4VProcess* p2 = aStep->GetPreStepPoint()->GetProcessDefinedStep();
00260   if (p)
00261     LogDebug("MuonSimDebug") <<" POST PROCESS = "<<p->GetProcessName()<<std::endl;
00262   if (p2)
00263     LogDebug("MuonSimDebug") <<" PRE  PROCESS = "<<p2->GetProcessName()<<std::endl;
00264   LogDebug("MuonSimDebug") << "newhit theta " << theThetaAtEntry<<std::endl;
00265   LogDebug("MuonSimDebug") << "newhit phi   " << thePhiAtEntry<<std::endl;
00266   LogDebug("MuonSimDebug") << "newhit pabs  " << thePabs<<std::endl;
00267   LogDebug("MuonSimDebug") << "newhit tof   " << theTof<<std::endl;
00268   LogDebug("MuonSimDebug") << "newhit track " << theTrackID<<std::endl;
00269   LogDebug("MuonSimDebug") << "newhit entry " << theEntryPoint<<std::endl;
00270   LogDebug("MuonSimDebug") << "newhit exit  " << theExitPoint<<std::endl;
00271   LogDebug("MuonSimDebug") << "newhit eloss " << theEnergyLoss << std::endl;
00272   LogDebug("MuonSimDebug") << "newhit detid " << theDetUnitId<<std::endl;
00273   LogDebug("MuonSimDebug") << "newhit delta " << (theExitPoint-theEntryPoint)<<std::endl;
00274   LogDebug("MuonSimDebug") << "newhit deltu " << (theExitPoint-theEntryPoint).unit();
00275   LogDebug("MuonSimDebug") << " " << (theExitPoint-theEntryPoint).mag()<<std::endl;
00276   LogDebug("MuonSimDebug") << "newhit glob  " << theGlobalEntry<<std::endl;
00277   LogDebug("MuonSimDebug") << "newhit dpos  " << theGlobalPos<<std::endl;
00278   LogDebug("MuonSimDebug") << "newhit drot  " << std::endl;
00279   //  theGlobalRot->print(LogDebug("MuonSimDebug"));
00280 
00281 
00282   //
00283   //----- SimTracks: Make it persistent?
00284   //
00285   int thePID = theTrack->GetDefinition()->GetPDGEncoding();
00286   LogDebug("MuonSimDebug") << " checking simtrack " << thePID << " " << thePabs << " STenergyPersistentCut " << STenergyPersistentCut << std::endl;
00287 
00288   if( thePabs*GeV > STenergyPersistentCut 
00289       || ( abs(thePID) == 13 && STallMuonsPersistent ) ){
00290     TrackInformation* info = getOrCreateTrackInformation(theTrack);
00291     LogDebug("MuonSimDebug") <<" track leaving hit in muons made selected for persistency"<<std::endl;
00292 
00293     info->storeTrack(true);
00294   }
00295      
00296 }
00297 
00298 void MuonSensitiveDetector::updateHit(G4Step * aStep){
00299   //  float thePabs             = aStep->GetPreStepPoint()->GetMomentum().mag()/GeV;
00300   //  Local3DPoint theEntryPoint= InitialStepPosition(aStep,LocalCoordinates);  
00301 
00302 
00303   Local3DPoint theExitPoint;
00304 
00305   if (detector->isBarrel()) {
00306     theExitPoint= toOrcaUnits(toOrcaRef(FinalStepPositionVsParent(aStep,1),aStep));  
00307   } else if (detector->isEndcap()) {
00308     // save local z at current level
00309     theExitPoint= toOrcaUnits(toOrcaRef(FinalStepPosition(aStep,LocalCoordinates),aStep));
00310     float zexit = theExitPoint.z();
00311     Local3DPoint tempExitPoint= toOrcaUnits(toOrcaRef(FinalStepPositionVsParent(aStep,4),aStep));
00312     theExitPoint  = Local3DPoint( tempExitPoint.x(), tempExitPoint.y(), zexit );
00313   } else {
00314     theExitPoint= toOrcaUnits(toOrcaRef(FinalStepPosition(aStep,LocalCoordinates),aStep)); 
00315   }
00316 
00317   float theEnergyLoss = aStep->GetTotalEnergyDeposit()/GeV;  
00318 
00319   if( theHit == 0 ){ 
00320     std::cerr << "!!ERRROR in MuonSensitiveDetector::updateHit. It is called when there is no hit " << std::endl;
00321   }
00322 
00323   theHit->updateExitPoint(theExitPoint);
00324   theHit->addEnergyLoss(theEnergyLoss);
00325 
00326   LogDebug("MuonSimDebug") <<"=== UPDATE ===============> ELOSS   = "<<theEnergyLoss<<" "
00327        <<thePV->GetLogicalVolume()->GetName()<<std::endl;
00328   const G4VProcess* p = aStep->GetPostStepPoint()->GetProcessDefinedStep();
00329   const G4VProcess* p2 = aStep->GetPreStepPoint()->GetProcessDefinedStep();
00330   if (p)
00331     LogDebug("MuonSimDebug") <<" POST PROCESS = "<<p->GetProcessName()<<std::endl;
00332   if (p2)
00333     LogDebug("MuonSimDebug") <<" PRE  PROCESS = "<<p2->GetProcessName()<<std::endl;
00334   LogDebug("MuonSimDebug") << "updhit exit  " << theExitPoint<<std::endl;
00335   LogDebug("MuonSimDebug") << "updhit eloss " << theHit->energyLoss() <<std::endl;
00336   LogDebug("MuonSimDebug") << "updhit detid " << theDetUnitId<<std::endl;
00337   LogDebug("MuonSimDebug") << "updhit delta " << (theExitPoint-theHit->entryPoint())<<std::endl;
00338   LogDebug("MuonSimDebug") << "updhit deltu " << (theExitPoint-theHit->entryPoint()).unit();
00339   LogDebug("MuonSimDebug") << " " << (theExitPoint-theHit->entryPoint()).mag()<<std::endl; 
00340 
00341 }
00342 
00343 void MuonSensitiveDetector::saveHit(){
00344 
00345   if (theHit) {
00346     if (printHits) {
00347       thePrinter->startNewSimHit(detector->name());
00348       thePrinter->printId(theHit->detUnitId());
00349       //      thePrinter->printTrack(theHit->trackId());
00350       //thePrinter->printPabs(theHit->pabs());
00351       //thePrinter->printEloss(theHit->energyLoss());
00352       thePrinter->printLocal(theHit->entryPoint(),theHit->exitPoint());
00353       thePrinter->printGlobal(theGlobalEntry);
00354     }
00355     slaveMuon->processHits(*theHit);
00356     // seems the hit does not want to be deleted
00357     // done by the hit collection?
00358     delete theHit;
00359     theHit = 0; //set it to 0, because you are checking that is 0
00360   }
00361 
00362 }
00363 
00364 TrackInformation* MuonSensitiveDetector::getOrCreateTrackInformation( const G4Track* gTrack)
00365 {
00366   G4VUserTrackInformation* temp = gTrack->GetUserInformation();
00367   if (temp == 0){
00368     std::cerr <<" ERROR: no G4VUserTrackInformation available"<<std::endl;
00369     abort();
00370   }else{
00371     TrackInformation* info = dynamic_cast<TrackInformation*>(temp);
00372     if (info ==0){
00373       std::cerr <<" ERROR: TkSimTrackSelection: the UserInformation does not appear to be a TrackInformation"<<std::endl;
00374       abort();
00375     }
00376     return info;
00377   }
00378 }
00379 
00380 void MuonSensitiveDetector::EndOfEvent(G4HCofThisEvent*)
00381 {
00382 //  TimeMe t("MuonSensitiveDetector::EndOfEvent", false);
00383  // LogDebug("MuonSimDebug") << "MuonSensitiveDetector::EndOfEvent saving last hit en event " << std::endl;
00384   saveHit();
00385 }
00386 
00387 
00388 void MuonSensitiveDetector::fillHits(edm::PSimHitContainer& c, std::string n){
00389   //
00390   // do it once for low, once for High
00391   //
00392 
00393   if (slaveMuon->name() == n) c=slaveMuon->hits();
00394 
00395 }
00396 
00397 std::vector<std::string> MuonSensitiveDetector::getNames(){
00398   std::vector<std::string> temp;
00399   temp.push_back(slaveMuon->name());
00400   return temp;
00401 }
00402 
00403 Local3DPoint MuonSensitiveDetector::InitialStepPositionVsParent(G4Step * currentStep, G4int levelsUp) {
00404   
00405   G4StepPoint * preStepPoint = currentStep->GetPreStepPoint();
00406   G4ThreeVector globalCoordinates = preStepPoint->GetPosition();
00407   
00408   G4TouchableHistory * theTouchable=(G4TouchableHistory *)
00409     (preStepPoint->GetTouchable());
00410 
00411   G4int depth = theTouchable->GetHistory()->GetDepth();
00412   G4ThreeVector localCoordinates = theTouchable->GetHistory()
00413     ->GetTransform(depth-levelsUp).TransformPoint(globalCoordinates);
00414   
00415   return ConvertToLocal3DPoint(localCoordinates); 
00416 }
00417  
00418 Local3DPoint MuonSensitiveDetector::FinalStepPositionVsParent(G4Step * currentStep, G4int levelsUp) {
00419   
00420   G4StepPoint * postStepPoint = currentStep->GetPostStepPoint();
00421   G4StepPoint * preStepPoint  = currentStep->GetPreStepPoint();
00422   G4ThreeVector globalCoordinates = postStepPoint->GetPosition();
00423     
00424   G4TouchableHistory * theTouchable = (G4TouchableHistory *)
00425     (preStepPoint->GetTouchable());
00426 
00427   G4int depth = theTouchable->GetHistory()->GetDepth();
00428   G4ThreeVector localCoordinates = theTouchable->GetHistory()
00429     ->GetTransform(depth-levelsUp).TransformPoint(globalCoordinates);
00430 
00431   return ConvertToLocal3DPoint(localCoordinates); 
00432 }