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 " << HCE;
41  theHC = new FiberG4HitsCollection(GetName(), collectionName[0]);
42  if (theHCID < 0)
43  theHCID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
44  HCE->AddHitsCollection(theHCID, theHC);
45  edm::LogVerbatim("FiberSim") << "FiberSD : Add hit collectrion for " << collectionName[0] << ":" << theHCID << ":"
46  << theHC;
47 }
48 
49 G4bool FiberSD::ProcessHits(G4Step* aStep, G4TouchableHistory*) {
50  //std::vector<HFShower::Hit> hits = theShower->getHits(aStep);
51  double zoffset = 1000;
52  std::vector<HFShower::Hit> hits = theShower->getHits(aStep, true, zoffset);
53 
54  if (!hits.empty()) {
55  std::vector<HFShowerPhoton> thePE;
56  for (unsigned int i = 0; i < hits.size(); i++) {
57  //edm::LogVerbatim("FiberSim") << "FiberSD :hit position z " << hits[i].position.z();
59  hits[i].position.x(), hits[i].position.y(), hits[i].position.z(), hits[i].wavelength, hits[i].time);
60  thePE.push_back(pe);
61  }
62  int trackID = aStep->GetTrack()->GetTrackID();
63  G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
64  const G4VTouchable* touch = preStepPoint->GetTouchable();
65  G4LogicalVolume* lv = touch->GetVolume(0)->GetLogicalVolume();
66  int depth = (touch->GetReplicaNumber(0)) % 10;
67  int detID = setDetUnitId(aStep);
68  math::XYZPoint theHitPos(
69  preStepPoint->GetPosition().x(), preStepPoint->GetPosition().y(), preStepPoint->GetPosition().z());
70  //edm::LogVerbatim("FiberSim") << "FiberSD :presteppoint position z " << preStepPoint->GetPosition().z();
71 
72  FiberG4Hit* aHit = new FiberG4Hit(lv, detID, depth, trackID);
73 #ifdef EDM_ML_DEBUG
74  edm::LogVerbatim("FiberSim") << "FiberSD :hit size " << hits.size() << " npe" << aHit->npe();
75  edm::LogVerbatim("FiberSim") << "FiberSD :pre hit position " << aHit->hitPos();
76 #endif
77  aHit->setNpe(hits.size());
78  aHit->setPos(theHitPos);
79  aHit->setTime(preStepPoint->GetGlobalTime());
80  aHit->setPhoton(thePE);
81 #ifdef EDM_ML_DEBUG
82  edm::LogVerbatim("FiberSim") << "FiberSD :ShowerPhoton position " << thePE[0].x() << " " << thePE[0].y() << " "
83  << thePE[0].z();
84 
85  edm::LogVerbatim("FiberSim") << "FiberSD: Hit created at " << lv->GetName() << " DetID: " << aHit->towerId()
86  << " Depth: " << aHit->depth() << " Track ID: " << aHit->trackId()
87  << " Nb. of Cerenkov Photons: " << aHit->npe() << " Time: " << aHit->time() << " at "
88  << aHit->hitPos();
89  for (unsigned int i = 0; i < thePE.size(); i++)
90  edm::LogVerbatim("FiberSim") << "FiberSD: PE[" << i << "] " << thePE[i];
91 #endif
92  theHC->insert(aHit);
93  }
94  return true;
95 }
96 
97 void FiberSD::EndOfEvent(G4HCofThisEvent* HCE) {
98  edm::LogVerbatim("FiberSim") << "FiberSD: finds " << theHC->entries() << " hits";
99  clear();
100  edm::LogVerbatim("FiberSim") << "theHC entries = " << theHC->entries();
101 }
102 
103 void FiberSD::clear() {}
104 
106 
108 
109 void FiberSD::update(const BeginOfJob* job) {}
110 
112 
114 
115 void FiberSD::update(const ::EndOfEvent*) {}
116 
118 
119 uint32_t FiberSD::setDetUnitId(const G4Step* aStep) {
120  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
121  int fibre = (touch->GetReplicaNumber(1)) % 10;
122  int cell = (touch->GetReplicaNumber(2));
123  int tower = (touch->GetReplicaNumber(3));
124  return ((tower * 1000 + cell) * 10 + fibre);
125 }
126 
Log< level::Info, true > LogVerbatim
uint32_t setDetUnitId(const G4Step *) override
Definition: FiberSD.cc:119
void clearHits() override
Definition: FiberSD.cc:117
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:40
const HcalSimulationParameters * hcalsimpar() const
void update(const BeginOfJob *) override
This routine will be called when the appropriate signal arrives.
Definition: FiberSD.cc:109
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:49
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:127
static int position[264][3]
Definition: ReadPGInfo.cc:289
void DrawAll() override
Definition: FiberSD.cc:105
void EndOfEvent(G4HCofThisEvent *HCE) override
Definition: FiberSD.cc:97
void clear() override
Definition: FiberSD.cc:103
G4int theHCID
Definition: FiberSD.h:63
void PrintAll() override
Definition: FiberSD.cc:107