CMS 3D CMS Logo

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