CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_1/src/SimG4CMS/ShowerLibraryProducer/src/HFWedgeSD.cc

Go to the documentation of this file.
00001 #include "SimG4CMS/ShowerLibraryProducer/interface/HFWedgeSD.h"
00002 #include "DataFormats/Math/interface/Point3D.h"
00003 
00004 #include "G4VPhysicalVolume.hh"
00005 #include "G4PVPlacement.hh"
00006 #include "G4HCofThisEvent.hh"
00007 #include "G4TouchableHistory.hh"
00008 #include "G4Track.hh"
00009 #include "G4Step.hh"
00010 #include "G4VSolid.hh"
00011 #include "G4DynamicParticle.hh"
00012 #include "G4ParticleDefinition.hh"
00013 #include "G4SDManager.hh"
00014 #include "G4ios.hh"
00015 
00016 HFWedgeSD::HFWedgeSD(std::string name, const DDCompactView & cpv, 
00017                  SensitiveDetectorCatalog & clg, edm::ParameterSet const & p, 
00018                  const SimTrackManager* manager) :
00019   SensitiveCaloDetector(name, cpv, clg, p), theName(name),
00020   m_trackManager(manager), hcID(-1), theHC(0), currentHit(0) {
00021 
00022   collectionName.insert(name);
00023   LogDebug("FiberSim") << "***************************************************"
00024                        << "\n"
00025                        << "*                                                 *"
00026                        << "\n"
00027                        << "* Constructing a HFWedgeSD  with name " << GetName()
00028                        << "\n"
00029                        << "*                                                 *"
00030                        << "\n"
00031                        << "***************************************************";
00032   //
00033   // Now attach the right detectors (LogicalVolumes) to me
00034   //
00035   std::vector<std::string> lvNames = clg.logicalNames(name);
00036   this->Register();
00037   for (std::vector<std::string>::iterator it=lvNames.begin();
00038        it !=lvNames.end(); it++){
00039     this->AssignSD(*it);
00040     LogDebug("FiberSim") << "HFWedgeSD : Assigns SD to LV " << (*it);
00041   }
00042 }
00043 
00044 HFWedgeSD::~HFWedgeSD() {
00045   if (theHC)    delete theHC;
00046 }
00047 
00048 void HFWedgeSD::Initialize(G4HCofThisEvent * HCE) {
00049 
00050   LogDebug("FiberSim") << "HFWedgeSD : Initialize called for " << GetName();
00051   theHC = new HFShowerG4HitsCollection(GetName(), collectionName[0]);
00052   if (hcID<0)
00053     hcID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
00054   HCE->AddHitsCollection(hcID, theHC);
00055 
00056   clearHits();
00057 }
00058 
00059 G4bool HFWedgeSD::ProcessHits(G4Step * aStep, G4TouchableHistory*) {
00060 
00061   G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
00062   const G4VTouchable* touch = preStepPoint->GetTouchable();
00063   currentID = setDetUnitId(aStep);
00064   trackID   = aStep->GetTrack()->GetTrackID();
00065   edep      = aStep->GetTotalEnergyDeposit();
00066   time      = (preStepPoint->GetGlobalTime())/ns;
00067 
00068   globalPos = preStepPoint->GetPosition();
00069   localPos  = touch->GetHistory()->GetTopTransform().TransformPoint(globalPos);
00070   const G4DynamicParticle* particle =  aStep->GetTrack()->GetDynamicParticle();
00071   momDir    = particle->GetMomentumDirection();
00072 
00073   if (hitExists() == false && edep>0.) 
00074     currentHit = createNewHit();
00075 
00076   return true;
00077 }
00078 
00079 void HFWedgeSD::EndOfEvent(G4HCofThisEvent * HCE) {
00080  
00081   LogDebug("FiberSim") << "HFWedgeSD: Sees" << theHC->entries() << " hits";
00082   clear();
00083 }
00084 
00085 void HFWedgeSD::clear() {}
00086 
00087 void HFWedgeSD::DrawAll()  {}
00088 
00089 void HFWedgeSD::PrintAll() {}
00090 
00091 G4bool HFWedgeSD::hitExists() {
00092    
00093   // Update if in the same detector, time-slice and for same track   
00094   if (currentID == previousID) {
00095     updateHit(currentHit);
00096     return true;
00097   }
00098    
00099   std::map<int,HFShowerG4Hit*>::const_iterator it = hitMap.find(currentID);
00100   if (it != hitMap.end()) {
00101     updateHit(currentHit);
00102     return true;
00103   }
00104   
00105   return false;
00106 }
00107 
00108 HFShowerG4Hit* HFWedgeSD::createNewHit() {
00109 
00110   LogDebug("FiberSim") << "HFWedgeSD::CreateNewHit for ID " << currentID
00111                        << " Track " << trackID << " Edep: " << edep/MeV 
00112                        << " MeV; Time: " << time << " ns; Position (local) " 
00113                        << localPos << " (global ) " << globalPos 
00114                        << " direction " << momDir;
00115   HFShowerG4Hit* aHit = new HFShowerG4Hit;
00116   aHit->setHitId(currentID);
00117   aHit->setTrackId(trackID);
00118   aHit->setTime(time);
00119   aHit->setLocalPos(localPos);
00120   aHit->setGlobalPos(globalPos);
00121   aHit->setPrimMomDir(momDir);
00122   updateHit(aHit);
00123 
00124   theHC->insert(aHit);
00125   hitMap.insert(std::pair<int,HFShowerG4Hit*>(previousID,aHit));
00126 
00127   return aHit;
00128 }        
00129 
00130 void HFWedgeSD::updateHit(HFShowerG4Hit* aHit) {
00131 
00132   if (edep != 0) {
00133     aHit->updateEnergy(edep);
00134     LogDebug("FiberSim") << "HFWedgeSD: Add energy deposit in " << currentID 
00135                          << " edep " << edep/MeV << " MeV"; 
00136   }
00137   previousID = currentID;
00138 }
00139 
00140 void HFWedgeSD::clearHits() {
00141   hitMap.erase (hitMap.begin(), hitMap.end());
00142   previousID = -1;
00143 }
00144 
00145 uint32_t HFWedgeSD::setDetUnitId(G4Step* aStep) {
00146   const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
00147   return (touch->GetReplicaNumber(0));
00148 }
00149 
00150 void HFWedgeSD::fillHits(edm::PCaloHitContainer&, std::string) {}