File indexing completed on 2024-04-06 12:30:12
0001 #include "SimG4CMS/ShowerLibraryProducer/interface/FiberSD.h"
0002 #include "SimDataFormats/CaloHit/interface/HFShowerPhoton.h"
0003 #include "DataFormats/Math/interface/Point3D.h"
0004 #include "Geometry/HcalCommonData/interface/HcalDDDSimConstants.h"
0005 #include "Geometry/HcalCommonData/interface/HcalSimulationConstants.h"
0006 #include "Geometry/Records/interface/HcalSimNumberingRecord.h"
0007
0008 #include "G4VPhysicalVolume.hh"
0009 #include "G4PVPlacement.hh"
0010 #include "G4HCofThisEvent.hh"
0011 #include "G4TouchableHistory.hh"
0012 #include "G4Track.hh"
0013 #include "G4Step.hh"
0014 #include "G4VSolid.hh"
0015 #include "G4DynamicParticle.hh"
0016 #include "G4ParticleDefinition.hh"
0017 #include "G4SDManager.hh"
0018 #include "G4ios.hh"
0019
0020
0021
0022 FiberSD::FiberSD(const std::string& iname,
0023 const HcalSimulationConstants* hsps,
0024 const HcalDDDSimConstants* hdc,
0025 const SensitiveDetectorCatalog& clg,
0026 edm::ParameterSet const& p,
0027 const SimTrackManager* manager)
0028 : SensitiveCaloDetector(iname, clg), theShower(nullptr), theHCID(-1), theHC(nullptr) {
0029 edm::LogVerbatim("FiberSim") << "FiberSD : Instantiating for " << iname;
0030
0031 theShower = new HFShower(iname, hdc, hsps->hcalsimpar(), p, 1);
0032 }
0033
0034 FiberSD::~FiberSD() {
0035 delete theShower;
0036 delete theHC;
0037 }
0038
0039 void FiberSD::Initialize(G4HCofThisEvent* HCE) {
0040 edm::LogVerbatim("FiberSim") << "FiberSD : Initialize called for " << GetName() << " in collection " << HCE;
0041 theHC = new FiberG4HitsCollection(GetName(), collectionName[0]);
0042 if (theHCID < 0)
0043 theHCID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
0044 HCE->AddHitsCollection(theHCID, theHC);
0045 edm::LogVerbatim("FiberSim") << "FiberSD : Add hit collectrion for " << collectionName[0] << ":" << theHCID << ":"
0046 << theHC;
0047 }
0048
0049 G4bool FiberSD::ProcessHits(G4Step* aStep, G4TouchableHistory*) {
0050
0051 double zoffset = 1000;
0052 std::vector<HFShower::Hit> hits = theShower->getHits(aStep, true, zoffset);
0053
0054 if (!hits.empty()) {
0055 std::vector<HFShowerPhoton> thePE;
0056 for (unsigned int i = 0; i < hits.size(); i++) {
0057
0058 HFShowerPhoton pe = HFShowerPhoton(
0059 hits[i].position.x(), hits[i].position.y(), hits[i].position.z(), hits[i].wavelength, hits[i].time);
0060 thePE.push_back(pe);
0061 }
0062 int trackID = aStep->GetTrack()->GetTrackID();
0063 G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
0064 const G4VTouchable* touch = preStepPoint->GetTouchable();
0065 G4LogicalVolume* lv = touch->GetVolume(0)->GetLogicalVolume();
0066 int depth = (touch->GetReplicaNumber(0)) % 10;
0067 int detID = setDetUnitId(aStep);
0068 math::XYZPoint theHitPos(
0069 preStepPoint->GetPosition().x(), preStepPoint->GetPosition().y(), preStepPoint->GetPosition().z());
0070
0071
0072 FiberG4Hit* aHit = new FiberG4Hit(lv, detID, depth, trackID);
0073 #ifdef EDM_ML_DEBUG
0074 edm::LogVerbatim("FiberSim") << "FiberSD :hit size " << hits.size() << " npe" << aHit->npe();
0075 edm::LogVerbatim("FiberSim") << "FiberSD :pre hit position " << aHit->hitPos();
0076 #endif
0077 aHit->setNpe(hits.size());
0078 aHit->setPos(theHitPos);
0079 aHit->setTime(preStepPoint->GetGlobalTime());
0080 aHit->setPhoton(thePE);
0081 #ifdef EDM_ML_DEBUG
0082 edm::LogVerbatim("FiberSim") << "FiberSD :ShowerPhoton position " << thePE[0].x() << " " << thePE[0].y() << " "
0083 << thePE[0].z();
0084
0085 edm::LogVerbatim("FiberSim") << "FiberSD: Hit created at " << lv->GetName() << " DetID: " << aHit->towerId()
0086 << " Depth: " << aHit->depth() << " Track ID: " << aHit->trackId()
0087 << " Nb. of Cerenkov Photons: " << aHit->npe() << " Time: " << aHit->time() << " at "
0088 << aHit->hitPos();
0089 for (unsigned int i = 0; i < thePE.size(); i++)
0090 edm::LogVerbatim("FiberSim") << "FiberSD: PE[" << i << "] " << thePE[i];
0091 #endif
0092 theHC->insert(aHit);
0093 }
0094 return true;
0095 }
0096
0097 void FiberSD::EndOfEvent(G4HCofThisEvent* HCE) {
0098 edm::LogVerbatim("FiberSim") << "FiberSD: finds " << theHC->entries() << " hits";
0099 clear();
0100 edm::LogVerbatim("FiberSim") << "theHC entries = " << theHC->entries();
0101 }
0102
0103 void FiberSD::clear() {}
0104
0105 void FiberSD::DrawAll() {}
0106
0107 void FiberSD::PrintAll() {}
0108
0109 void FiberSD::update(const BeginOfJob* job) {}
0110
0111 void FiberSD::update(const BeginOfRun*) {}
0112
0113 void FiberSD::update(const BeginOfEvent*) {}
0114
0115 void FiberSD::update(const ::EndOfEvent*) {}
0116
0117 void FiberSD::clearHits() {}
0118
0119 uint32_t FiberSD::setDetUnitId(const G4Step* aStep) {
0120 const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
0121 int fibre = (touch->GetReplicaNumber(1)) % 10;
0122 int cell = (touch->GetReplicaNumber(2));
0123 int tower = (touch->GetReplicaNumber(3));
0124 return ((tower * 1000 + cell) * 10 + fibre);
0125 }
0126
0127 void FiberSD::fillHits(edm::PCaloHitContainer&, const std::string&) {}