CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch13/src/SimG4CMS/Forward/src/BscSD.cc

Go to the documentation of this file.
00001 
00002 // File: BscSD.cc
00003 // Date: 02.2006
00004 // Description: Sensitive Detector class for Bsc
00005 // Modifications: 
00007  
00008 //#include "Geometry/Vector/interface/LocalPoint.h"
00009 //#include "Geometry/Vector/interface/LocalVector.h"
00010 
00011 #include "SimG4Core/Notification/interface/TrackInformation.h"
00012 #include "SimG4Core/Notification/interface/G4TrackToParticleID.h"
00013 #include "SimG4Core/Physics/interface/G4ProcessTypeEnumerator.h"
00014 
00015 #include "SimDataFormats/SimHitMaker/interface/TrackingSlaveSD.h"
00016 #include "SimDataFormats/TrackingHit/interface/UpdatablePSimHit.h"
00017 
00018 
00019 #include "SimG4CMS/Forward/interface/BscSD.h"
00020 #include "SimG4CMS/Forward/interface/BscG4Hit.h"
00021 #include "SimG4CMS/Forward/interface/BscG4HitCollection.h"
00022 #include "SimG4CMS/Forward/interface/BscNumberingScheme.h"
00023 
00024 #include "FWCore/Framework/interface/ESHandle.h"
00025 #include "FWCore/Framework/interface/EventSetup.h"
00026 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00027 
00028 #include "SimG4Core/Notification/interface/TrackInformation.h"
00029 
00030 #include "G4Track.hh"
00031 #include "G4SDManager.hh"
00032 #include "G4VProcess.hh"
00033 #include "G4EventManager.hh"
00034 #include "G4Step.hh"
00035 #include "G4ParticleTable.hh"
00036 
00037 #include <string>
00038 #include <vector>
00039 #include <iostream>
00040 
00041 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00042 
00043 #define debug
00044 //-------------------------------------------------------------------
00045 BscSD::BscSD(std::string name, const DDCompactView & cpv,
00046              SensitiveDetectorCatalog & clg, 
00047              edm::ParameterSet const & p, const SimTrackManager* manager) :
00048   SensitiveTkDetector(name, cpv, clg, p), numberingScheme(0), name(name),
00049   hcID(-1), theHC(0), theManager(manager), currentHit(0), theTrack(0), 
00050   currentPV(0), unitID(0),  previousUnitID(0), preStepPoint(0), 
00051   postStepPoint(0), eventno(0){
00052     
00053   //Add Bsc Sentitive Detector Name
00054   collectionName.insert(name);
00055     
00056     
00057   //Parameters
00058   edm::ParameterSet m_p = p.getParameter<edm::ParameterSet>("BscSD");
00059   int verbn = m_p.getUntrackedParameter<int>("Verbosity");
00060   //int verbn = 1;
00061     
00062   SetVerboseLevel(verbn);
00063   LogDebug("BscSim") 
00064     << "*******************************************************\n"
00065     << "*                                                     *\n"
00066     << "* Constructing a BscSD  with name " << name << "\n"
00067     << "*                                                     *\n"
00068     << "*******************************************************";
00069     
00070     
00071   slave  = new TrackingSlaveSD(name);
00072     
00073   //
00074   // attach detectors (LogicalVolumes)
00075   //
00076   std::vector<std::string> lvNames = clg.logicalNames(name);
00077 
00078   this->Register();
00079 
00080   for (std::vector<std::string>::iterator it=lvNames.begin();  
00081        it !=lvNames.end(); it++) {
00082     this->AssignSD(*it);
00083     edm::LogInfo("BscSim") << "BscSD : Assigns SD to LV " << (*it);
00084   }
00085     
00086   if      (name == "BSCHits") {
00087     if (verbn > 0) {
00088       edm::LogInfo("BscSim") << "name = BSCHits and  new BscNumberingSchem";
00089     }
00090     numberingScheme = new BscNumberingScheme() ;
00091   } else {
00092     edm::LogWarning("BscSim") << "BscSD: ReadoutName "<<name<<" not supported";
00093   }
00094   
00095   edm::LogInfo("BscSim") << "BscSD: Instantiation completed";
00096 }
00097 
00098 
00099 
00100 
00101 BscSD::~BscSD() { 
00102   //AZ:
00103   if (slave) delete slave; 
00104 
00105   if (numberingScheme)
00106     delete numberingScheme;
00107 
00108 }
00109 
00110 double BscSD::getEnergyDeposit(G4Step* aStep) {
00111   return aStep->GetTotalEnergyDeposit();
00112 }
00113 
00114 void BscSD::Initialize(G4HCofThisEvent * HCE) { 
00115 #ifdef debug
00116   LogDebug("BscSim") << "BscSD : Initialize called for " << name << std::endl;
00117 #endif
00118 
00119   theHC = new BscG4HitCollection(name, collectionName[0]);
00120   if (hcID<0) 
00121     hcID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
00122   HCE->AddHitsCollection(hcID, theHC);
00123 
00124   tsID   = -2;
00125   primID = -2;
00126 
00128 }
00129 
00130 
00131 bool BscSD::ProcessHits(G4Step * aStep, G4TouchableHistory * ) {
00132 
00133   if (aStep == NULL) {
00134     return true;
00135   } else {
00136     GetStepInfo(aStep);
00137     //   LogDebug("BscSim") << edeposit <<std::endl;
00138 
00139     //AZ
00140 #ifdef debug
00141     LogDebug("BscSim") << "BscSD :  number of hits = " << theHC->entries() << std::endl;
00142 #endif
00143 
00144     if (HitExists() == false && edeposit>0. ){ 
00145       CreateNewHit();
00146       return true;
00147     }
00148   }
00149   return true;
00150 } 
00151 
00152 void BscSD::GetStepInfo(G4Step* aStep) {
00153   
00154   preStepPoint = aStep->GetPreStepPoint(); 
00155   postStepPoint= aStep->GetPostStepPoint(); 
00156   theTrack     = aStep->GetTrack();   
00157   hitPoint     = preStepPoint->GetPosition();   
00158   currentPV    = preStepPoint->GetPhysicalVolume();
00159   hitPointExit = postStepPoint->GetPosition();  
00160 
00161   hitPointLocal = preStepPoint->GetTouchable()->GetHistory()->GetTopTransform().TransformPoint(hitPoint);
00162   hitPointLocalExit = preStepPoint->GetTouchable()->GetHistory()->GetTopTransform().TransformPoint(hitPointExit);
00163 
00164 
00165   G4int particleCode = theTrack->GetDefinition()->GetPDGEncoding();
00166   LogDebug("BscSim") <<  "  BscSD :particleType =  " << theTrack->GetDefinition()->GetParticleName() <<std::endl;
00167   if (particleCode == emPDG ||
00168       particleCode == epPDG ||
00169       particleCode == gammaPDG ) {
00170     edepositEM  = getEnergyDeposit(aStep); edepositHAD = 0.;
00171   } else {
00172     edepositEM  = 0.; edepositHAD = getEnergyDeposit(aStep);
00173   }
00174   edeposit = aStep->GetTotalEnergyDeposit();
00175   tSlice    = (postStepPoint->GetGlobalTime() )/nanosecond;
00176   tSliceID  = (int) tSlice;
00177   unitID    = setDetUnitId(aStep);
00178 #ifdef debug
00179   LogDebug("BscSim") << "unitID=" << unitID <<std::endl;
00180 #endif
00181   primaryID    = theTrack->GetTrackID();
00182   //  Position     = hitPoint;
00183   Pabs         = aStep->GetPreStepPoint()->GetMomentum().mag()/GeV;
00184   Tof          = aStep->GetPostStepPoint()->GetGlobalTime()/nanosecond;
00185   Eloss        = aStep->GetTotalEnergyDeposit()/GeV;
00186   ParticleType = theTrack->GetDefinition()->GetPDGEncoding();      
00187   ThetaAtEntry = aStep->GetPreStepPoint()->GetPosition().theta()/deg;
00188   PhiAtEntry   = aStep->GetPreStepPoint()->GetPosition().phi()/deg;
00189 
00190   ParentId = theTrack->GetParentID();
00191   Vx = theTrack->GetVertexPosition().x();
00192   Vy = theTrack->GetVertexPosition().y();
00193   Vz = theTrack->GetVertexPosition().z();
00194   X  = hitPoint.x();
00195   Y  = hitPoint.y();
00196   Z  = hitPoint.z();
00197 }
00198 
00199 uint32_t BscSD::setDetUnitId(G4Step * aStep) { 
00200 
00201   return (numberingScheme == 0 ? 0 : numberingScheme->getUnitID(aStep));
00202 }
00203 
00204 
00205 G4bool BscSD::HitExists() {
00206   if (primaryID<1) {
00207     edm::LogWarning("BscSim") << "***** BscSD error: primaryID = " 
00208                                   << primaryID
00209                                   << " maybe detector name changed";
00210   }
00211 
00212   // Update if in the same detector, time-slice and for same track   
00213   //  if (primaryID == primID && tSliceID == tsID && unitID==previousUnitID) {
00214   if (tSliceID == tsID && unitID==previousUnitID) {
00215     //AZ:
00216     UpdateHit();
00217     return true;
00218   }
00219   // Reset entry point for new primary
00220   if (primaryID != primID)
00221     ResetForNewPrimary();
00222    
00223   //look in the HitContainer whether a hit with the same primID, unitID,
00224   //tSliceID already exists:
00225    
00226   G4bool found = false;
00227 
00228   //    LogDebug("BscSim") << "BscSD: HCollection=  " << theHC->entries()    <<std::endl;
00229   
00230   for (int j=0; j<theHC->entries()&&!found; j++) {
00231     BscG4Hit* aPreviousHit = (*theHC)[j];
00232     if (aPreviousHit->getTrackID()     == primaryID &&
00233         aPreviousHit->getTimeSliceID() == tSliceID  &&
00234         aPreviousHit->getUnitID()      == unitID       ) {
00235       //AZ:
00236       currentHit = aPreviousHit;
00237       found      = true;
00238     }
00239   }          
00240 
00241   if (found) {
00242     //AZ:
00243     UpdateHit();
00244     return true;
00245   } else {
00246     return false;
00247   }    
00248 }
00249 
00250 
00251 void BscSD::ResetForNewPrimary() {
00252   
00253   entrancePoint  = SetToLocal(hitPoint);
00254   exitPoint      = SetToLocalExit(hitPointExit);
00255   incidentEnergy = preStepPoint->GetKineticEnergy();
00256 
00257 }
00258 
00259 
00260 void BscSD::StoreHit(BscG4Hit* hit){
00261 
00262   if (primID<0) return;
00263   if (hit == 0 ) {
00264     edm::LogWarning("BscSim") << "BscSD: hit to be stored is NULL !!";
00265     return;
00266   }
00267 
00268   theHC->insert( hit );
00269 }
00270 
00271 
00272 void BscSD::CreateNewHit() {
00273 
00274 #ifdef debug
00275   LogDebug("BscSim") << "BscSD CreateNewHit for"
00276                      << " PV "     << currentPV->GetName()
00277                      << " PVid = " << currentPV->GetCopyNo()
00278                      << " Unit "   << unitID <<std::endl;
00279   LogDebug("BscSim") << " primary "    << primaryID
00280                      << " time slice " << tSliceID 
00281                      << " For Track  " << theTrack->GetTrackID()
00282                      << " which is a " << theTrack->GetDefinition()->GetParticleName();
00283            
00284   if (theTrack->GetTrackID()==1) {
00285     LogDebug("BscSim") << " of energy "     << theTrack->GetTotalEnergy();
00286   } else {
00287     LogDebug("BscSim") << " daughter of part. " << theTrack->GetParentID();
00288   }
00289 
00290   LogDebug("BscSim")  << " and created by " ;
00291   if (theTrack->GetCreatorProcess()!=NULL)
00292     LogDebug("BscSim") << theTrack->GetCreatorProcess()->GetProcessName() ;
00293   else 
00294     LogDebug("BscSim") << "NO process";
00295   LogDebug("BscSim") << std::endl;
00296 #endif          
00297     
00298   currentHit = new BscG4Hit;
00299   currentHit->setTrackID(primaryID);
00300   currentHit->setTimeSlice(tSlice);
00301   currentHit->setUnitID(unitID);
00302   currentHit->setIncidentEnergy(incidentEnergy);
00303 
00304   currentHit->setPabs(Pabs);
00305   currentHit->setTof(Tof);
00306   currentHit->setEnergyLoss(Eloss);
00307   currentHit->setParticleType(ParticleType);
00308   currentHit->setThetaAtEntry(ThetaAtEntry);
00309   currentHit->setPhiAtEntry(PhiAtEntry);
00310 
00311   currentHit->setEntry(hitPoint);
00312 
00313   currentHit->setEntryLocalP(hitPointLocal);
00314   currentHit->setExitLocalP(hitPointLocalExit);
00315 
00316   currentHit->setParentId(ParentId);
00317   currentHit->setVx(Vx);
00318   currentHit->setVy(Vy);
00319   currentHit->setVz(Vz);
00320 
00321   currentHit->setX(X);
00322   currentHit->setY(Y);
00323   currentHit->setZ(Z);
00324 
00325   UpdateHit();
00326   
00327   StoreHit(currentHit);
00328 }        
00329  
00330 
00331 void BscSD::UpdateHit() {
00332 
00333   if (Eloss > 0.) {
00334     currentHit->addEnergyDeposit(edepositEM,edepositHAD);
00335 
00336 #ifdef debug
00337     LogDebug("BscSim") << "updateHit: add eloss " << Eloss <<std::endl;
00338     LogDebug("BscSim") << "CurrentHit=" << currentHit
00339                        << ", PostStepPoint=" << postStepPoint->GetPosition();
00340 #endif
00341     //AZ
00342     currentHit->setEnergyLoss(Eloss);
00343   }  
00344 
00345   // buffer for next steps:
00346   tsID           = tSliceID;
00347   primID         = primaryID;
00348   previousUnitID = unitID;
00349 }
00350 
00351 
00352 G4ThreeVector BscSD::SetToLocal(G4ThreeVector global){
00353 
00354   const G4VTouchable* touch= preStepPoint->GetTouchable();
00355   theEntryPoint = touch->GetHistory()->GetTopTransform().TransformPoint(global);
00356   return theEntryPoint;  
00357 }
00358      
00359 
00360 G4ThreeVector BscSD::SetToLocalExit(G4ThreeVector globalPoint){
00361 
00362   const G4VTouchable* touch= postStepPoint->GetTouchable();
00363   theExitPoint = touch->GetHistory()->GetTopTransform().TransformPoint(globalPoint);
00364   return theExitPoint;  
00365 }
00366      
00367 
00368 void BscSD::EndOfEvent(G4HCofThisEvent* ) {
00369 
00370   // here we loop over transient hits and make them persistent
00371   for (int j=0; j<theHC->entries(); j++) {
00372     //AZ:
00373     BscG4Hit* aHit = (*theHC)[j];
00374     LogDebug("BscSim") << "hit number" << j << "unit ID = "<<aHit->getUnitID()<< "\n";
00375     LogDebug("BscSim") << "entry z " << aHit->getEntry().z()<< "\n";
00376     LogDebug("BscSim") << "entr theta " << aHit->getThetaAtEntry()<< "\n";
00377 
00378     Local3DPoint locExitPoint(0,0,0);
00379     Local3DPoint locEntryPoint(aHit->getEntry().x(),
00380                                aHit->getEntry().y(),
00381                                aHit->getEntry().z());
00382     slave->processHits(PSimHit(locEntryPoint,locExitPoint,
00383                                aHit->getPabs(),
00384                                aHit->getTof(),
00385                                aHit->getEnergyLoss(),
00386                                aHit->getParticleType(),
00387                                aHit->getUnitID(),
00388                                aHit->getTrackID(),
00389                                aHit->getThetaAtEntry(),
00390                                aHit->getPhiAtEntry()));
00391   }
00392   Summarize();
00393 }
00394      
00395 
00396 void BscSD::Summarize() {
00397 }
00398 
00399 
00400 void BscSD::clear() {
00401 } 
00402 
00403 
00404 void BscSD::DrawAll() {
00405 } 
00406 
00407 
00408 void BscSD::PrintAll() {
00409   LogDebug("BscSim") << "BscSD: Collection " << theHC->GetName() << "\n";
00410   theHC->PrintAllHits();
00411 } 
00412 
00413 
00414 void BscSD::fillHits(edm::PSimHitContainer& c, std::string n) {
00415   if (slave->name() == n) c=slave->hits();
00416 }
00417 
00418 void BscSD::update (const BeginOfEvent * i) {
00419   LogDebug("BscSim") << " Dispatched BeginOfEvent for " << GetName()
00420                        << " !" ;
00421    clearHits();
00422    eventno = (*i)()->GetEventID();
00423 }
00424 
00425 void BscSD::update(const BeginOfRun *) {
00426 
00427   G4ParticleTable * theParticleTable = G4ParticleTable::GetParticleTable();
00428   G4String particleName;
00429   emPDG = theParticleTable->FindParticle(particleName="e-")->GetPDGEncoding();
00430   epPDG = theParticleTable->FindParticle(particleName="e+")->GetPDGEncoding();
00431   gammaPDG = theParticleTable->FindParticle(particleName="gamma")->GetPDGEncoding();
00432 
00433 } 
00434 
00435 void BscSD::update (const ::EndOfEvent*) {
00436 }
00437 
00438 void BscSD::clearHits(){
00439   //AZ:
00440   slave->Initialize();
00441 }
00442 
00443 std::vector<std::string> BscSD::getNames(){
00444   std::vector<std::string> temp;
00445   temp.push_back(slave->name());
00446   return temp;
00447 }
00448