File indexing completed on 2024-05-10 02:21:22
0001 #include "SimG4CMS/ShowerLibraryProducer/interface/HFWedgeSD.h"
0002 #include "DataFormats/Math/interface/Point3D.h"
0003
0004 #include "G4VPhysicalVolume.hh"
0005 #include "G4PVPlacement.hh"
0006 #include "G4HCofThisEvent.hh"
0007 #include "G4TouchableHistory.hh"
0008 #include "G4Track.hh"
0009 #include "G4Step.hh"
0010 #include "G4VSolid.hh"
0011 #include "G4DynamicParticle.hh"
0012 #include "G4ParticleDefinition.hh"
0013 #include "G4SDManager.hh"
0014 #include "G4ios.hh"
0015
0016 #include "G4PhysicalConstants.hh"
0017
0018 #include <CLHEP/Units/SystemOfUnits.h>
0019 using CLHEP::ns;
0020
0021
0022
0023 HFWedgeSD::HFWedgeSD(const std::string& iname, const SensitiveDetectorCatalog& clg, const SimTrackManager* manager)
0024 : SensitiveCaloDetector(iname, clg), hcID(-1), theHC(nullptr), currentHit(nullptr) {
0025 edm::LogVerbatim("FiberSim") << "HFWedgeSD : Instantiated for " << iname;
0026 }
0027
0028 HFWedgeSD::~HFWedgeSD() { delete theHC; }
0029
0030 void HFWedgeSD::Initialize(G4HCofThisEvent* HCE) {
0031 edm::LogVerbatim("FiberSim") << "HFWedgeSD : Initialize called for " << GetName() << " in collection " << HCE;
0032 theHC = new HFShowerG4HitsCollection(GetName(), collectionName[0]);
0033 if (hcID < 0)
0034 hcID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
0035 HCE->AddHitsCollection(hcID, theHC);
0036 edm::LogVerbatim("FiberSim") << "HFWedgeSD : Add hit collectrion for " << collectionName[0] << ":" << hcID << ":"
0037 << theHC;
0038
0039 clearHits();
0040 }
0041
0042 G4bool HFWedgeSD::ProcessHits(G4Step* aStep, G4TouchableHistory*) {
0043 G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
0044 const G4VTouchable* touch = preStepPoint->GetTouchable();
0045 currentID = setDetUnitId(aStep);
0046 trackID = aStep->GetTrack()->GetTrackID();
0047 edep = aStep->GetTotalEnergyDeposit();
0048 time = (preStepPoint->GetGlobalTime()) / ns;
0049
0050 globalPos = preStepPoint->GetPosition();
0051 localPos = touch->GetHistory()->GetTopTransform().TransformPoint(globalPos);
0052 const G4DynamicParticle* particle = aStep->GetTrack()->GetDynamicParticle();
0053 momDir = particle->GetMomentumDirection();
0054
0055 if (hitExists() == false && edep > 0.)
0056 currentHit = createNewHit();
0057
0058 return true;
0059 }
0060
0061 void HFWedgeSD::EndOfEvent(G4HCofThisEvent* HCE) {
0062 edm::LogVerbatim("FiberSim") << "HFWedgeSD: Sees" << theHC->entries() << " hits";
0063 clear();
0064 }
0065
0066 void HFWedgeSD::clear() {}
0067
0068 void HFWedgeSD::DrawAll() {}
0069
0070 void HFWedgeSD::PrintAll() {}
0071
0072 G4bool HFWedgeSD::hitExists() {
0073
0074 if (currentID == previousID) {
0075 updateHit(currentHit);
0076 return true;
0077 }
0078
0079 std::map<int, HFShowerG4Hit*>::const_iterator it = hitMap.find(currentID);
0080 if (it != hitMap.end()) {
0081 updateHit(currentHit);
0082 return true;
0083 }
0084
0085 return false;
0086 }
0087
0088 HFShowerG4Hit* HFWedgeSD::createNewHit() {
0089 #ifdef EDM_ML_DEBUG
0090 edm::LogVerbatim("FiberSim") << "HFWedgeSD::CreateNewHit for ID " << currentID << " Track " << trackID
0091 << " Edep: " << edep / CLHEP::MeV << " MeV; Time: " << time << " ns; Position (local) "
0092 << localPos << " (global ) " << globalPos << " direction " << momDir;
0093 #endif
0094 HFShowerG4Hit* aHit = new HFShowerG4Hit;
0095 aHit->setHitId(currentID);
0096 aHit->setTrackId(trackID);
0097 aHit->setTime(time);
0098 aHit->setLocalPos(localPos);
0099 aHit->setGlobalPos(globalPos);
0100 aHit->setPrimMomDir(momDir);
0101 updateHit(aHit);
0102
0103 theHC->insert(aHit);
0104 hitMap.insert(std::pair<int, HFShowerG4Hit*>(previousID, aHit));
0105
0106 return aHit;
0107 }
0108
0109 void HFWedgeSD::updateHit(HFShowerG4Hit* aHit) {
0110 if (edep != 0) {
0111 aHit->updateEnergy(edep);
0112 #ifdef EDM_ML_DEBUG
0113 edm::LogVerbatim("FiberSim") << "HFWedgeSD: Add energy deposit in " << currentID << " edep " << edep / CLHEP::MeV
0114 << " MeV";
0115 #endif
0116 }
0117 previousID = currentID;
0118 }
0119
0120 void HFWedgeSD::clearHits() {
0121 hitMap.erase(hitMap.begin(), hitMap.end());
0122 previousID = -1;
0123 }
0124
0125 uint32_t HFWedgeSD::setDetUnitId(const G4Step* aStep) {
0126 const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
0127 return (touch->GetReplicaNumber(0));
0128 }
0129
0130 void HFWedgeSD::fillHits(edm::PCaloHitContainer&, const std::string&) {}