CMS 3D CMS Logo

FiberSD.cc
Go to the documentation of this file.
7 
8 #include "G4VPhysicalVolume.hh"
9 #include "G4PVPlacement.hh"
10 #include "G4HCofThisEvent.hh"
11 #include "G4TouchableHistory.hh"
12 #include "G4Track.hh"
13 #include "G4Step.hh"
14 #include "G4VSolid.hh"
15 #include "G4DynamicParticle.hh"
16 #include "G4ParticleDefinition.hh"
17 #include "G4SDManager.hh"
18 #include "G4ios.hh"
19 
20 //#define EDM_ML_DEBUG
21 
23  const HcalSimulationConstants* hsps,
24  const HcalDDDSimConstants* hdc,
25  const SensitiveDetectorCatalog& clg,
26  edm::ParameterSet const& p,
27  const SimTrackManager* manager)
28  : SensitiveCaloDetector(iname, clg), theShower(nullptr), theHCID(-1), theHC(nullptr) {
29  edm::LogVerbatim("FiberSim") << "FiberSD : Instantiating for " << iname;
30  // Get pointer to HcalDDDConstants and HcalSimulationConstants
31  theShower = new HFShower(iname, hdc, hsps->hcalsimpar(), p, 1);
32 }
33 
35  delete theShower;
36  delete theHC;
37 }
38 
39 void FiberSD::Initialize(G4HCofThisEvent* HCE) {
40  edm::LogVerbatim("FiberSim") << "FiberSD : Initialize called for " << GetName() << " in collection "
41  << HCE;
42  theHC = new FiberG4HitsCollection(GetName(), collectionName[0]);
43  if (theHCID < 0)
44  theHCID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
45  HCE->AddHitsCollection(theHCID, theHC);
46  edm::LogVerbatim("FiberSim") << "FiberSD : Add hit collectrion for " << collectionName[0] << ":"
47  << theHCID << ":" << theHC;
48 }
49 
50 G4bool FiberSD::ProcessHits(G4Step* aStep, G4TouchableHistory*) {
51  //std::vector<HFShower::Hit> hits = theShower->getHits(aStep);
52  double zoffset = 1000;
53  std::vector<HFShower::Hit> hits = theShower->getHits(aStep, true, zoffset);
54 
55  if (!hits.empty()) {
56  std::vector<HFShowerPhoton> thePE;
57  for (unsigned int i = 0; i < hits.size(); i++) {
58  //edm::LogVerbatim("FiberSim") << "FiberSD :hit position z " << hits[i].position.z();
60  hits[i].position.x(), hits[i].position.y(), hits[i].position.z(), hits[i].wavelength, hits[i].time);
61  thePE.push_back(pe);
62  }
63  int trackID = aStep->GetTrack()->GetTrackID();
64  G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
65  const G4VTouchable* touch = preStepPoint->GetTouchable();
66  G4LogicalVolume* lv = touch->GetVolume(0)->GetLogicalVolume();
67  int depth = (touch->GetReplicaNumber(0)) % 10;
68  int detID = setDetUnitId(aStep);
69  math::XYZPoint theHitPos(
70  preStepPoint->GetPosition().x(), preStepPoint->GetPosition().y(), preStepPoint->GetPosition().z());
71  //edm::LogVerbatim("FiberSim") << "FiberSD :presteppoint position z " << preStepPoint->GetPosition().z();
72 
73  FiberG4Hit* aHit = new FiberG4Hit(lv, detID, depth, trackID);
74 #ifdef EDM_ML_DEBUG
75  edm::LogVerbatim("FiberSim") << "FiberSD :hit size " << hits.size() << " npe" << aHit->npe();
76  edm::LogVerbatim("FiberSim") << "FiberSD :pre hit position " << aHit->hitPos();
77 #endif
78  aHit->setNpe(hits.size());
79  aHit->setPos(theHitPos);
80  aHit->setTime(preStepPoint->GetGlobalTime());
81  aHit->setPhoton(thePE);
82 #ifdef EDM_ML_DEBUG
83  edm::LogVerbatim("FiberSim") << "FiberSD :ShowerPhoton position " << thePE[0].x() << " "
84  << thePE[0].y() << " " << thePE[0].z();
85 
86  edm::LogVerbatim("FiberSim") << "FiberSD: Hit created at " << lv->GetName()
87  << " DetID: " << aHit->towerId() << " Depth: " << aHit->depth()
88  << " Track ID: " << aHit->trackId() << " Nb. of Cerenkov Photons: " << aHit->npe()
89  << " Time: " << aHit->time() << " at " << aHit->hitPos();
90  for (unsigned int i = 0; i < thePE.size(); i++)
91  edm::LogVerbatim("FiberSim") << "FiberSD: PE[" << i << "] " << thePE[i];
92 #endif
93  theHC->insert(aHit);
94  }
95  return true;
96 }
97 
98 void FiberSD::EndOfEvent(G4HCofThisEvent* HCE) {
99  edm::LogVerbatim("FiberSim") << "FiberSD: finds " << theHC->entries() << " hits";
100  clear();
101  edm::LogVerbatim("FiberSim") << "theHC entries = " << theHC->entries();
102 }
103 
104 void FiberSD::clear() {}
105 
107 
109 
110 void FiberSD::update(const BeginOfJob* job) {}
111 
113 
115 
116 void FiberSD::update(const ::EndOfEvent*) {}
117 
119 
120 uint32_t FiberSD::setDetUnitId(const G4Step* aStep) {
121  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
122  int fibre = (touch->GetReplicaNumber(1)) % 10;
123  int cell = (touch->GetReplicaNumber(2));
124  int tower = (touch->GetReplicaNumber(3));
125  return ((tower * 1000 + cell) * 10 + fibre);
126 }
127 
Log< level::Info, true > LogVerbatim
uint32_t setDetUnitId(const G4Step *) override
Definition: FiberSD.cc:120
void clearHits() override
Definition: FiberSD.cc:118
std::vector< PCaloHit > PCaloHitContainer
FiberSD(const std::string &, const HcalSimulationConstants *, const HcalDDDSimConstants *, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
Definition: FiberSD.cc:22
std::vector< Hit > getHits(const G4Step *aStep, double weight)
Definition: HFShower.cc:45
const HcalSimulationParameters * hcalsimpar() const
void update(const BeginOfJob *) override
This routine will be called when the appropriate signal arrives.
Definition: FiberSD.cc:110
void Initialize(G4HCofThisEvent *HCE) override
Definition: FiberSD.cc:39
HFShower * theShower
Definition: FiberSD.h:61
~FiberSD() override
Definition: FiberSD.cc:34
G4bool ProcessHits(G4Step *aStep, G4TouchableHistory *ROhist) override
Definition: FiberSD.cc:50
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
FiberG4HitsCollection * theHC
Definition: FiberSD.h:64
G4THitsCollection< FiberG4Hit > FiberG4HitsCollection
Definition: FiberG4Hit.h:54
void fillHits(edm::PCaloHitContainer &, const std::string &) override
Definition: FiberSD.cc:128
static int position[264][3]
Definition: ReadPGInfo.cc:289
void DrawAll() override
Definition: FiberSD.cc:106
void EndOfEvent(G4HCofThisEvent *HCE) override
Definition: FiberSD.cc:98
void clear() override
Definition: FiberSD.cc:104
G4int theHCID
Definition: FiberSD.h:63
void PrintAll() override
Definition: FiberSD.cc:108