CMS 3D CMS Logo

FiberSD.cc
Go to the documentation of this file.
9 
10 #include "G4VPhysicalVolume.hh"
11 #include "G4PVPlacement.hh"
12 #include "G4HCofThisEvent.hh"
13 #include "G4TouchableHistory.hh"
14 #include "G4Track.hh"
15 #include "G4Step.hh"
16 #include "G4VSolid.hh"
17 #include "G4DynamicParticle.hh"
18 #include "G4ParticleDefinition.hh"
19 #include "G4SDManager.hh"
20 #include "G4ios.hh"
21 
23  const edm::EventSetup& es,
24  const SensitiveDetectorCatalog& clg,
25  edm::ParameterSet const& p,
26  const SimTrackManager* manager)
27  : SensitiveCaloDetector(iname, es, clg, p),
28  m_trackManager(manager),
29  theShower(nullptr),
30  theHCID(-1),
31  theHC(nullptr) {
32  edm::LogVerbatim("FiberSim") << "FiberSD : Instantiating for " << iname;
33  // Get pointer to HcalDDDConstants and HcalSimulationConstants
35  es.get<HcalSimNumberingRecord>().get(hdsc);
36  if (!hdsc.isValid()) {
37  edm::LogError("FiberSim") << "FiberSD : Cannot find HcalSimulationConstants";
38  throw cms::Exception("Unknown", "FiberSD") << "Cannot find HcalSimulationConstants\n";
39  }
40  const HcalSimulationConstants* hsps = hdsc.product();
42  es.get<HcalSimNumberingRecord>().get(hdc);
43  if (hdc.isValid()) {
44  const HcalDDDSimConstants* hcalConstants = hdc.product();
45  theShower = new HFShower(iname, hcalConstants, hsps->hcalsimpar(), p, 1);
46  } else {
47  edm::LogError("FiberSim") << "FiberSD : Cannot find HcalDDDSimConstant";
48  throw cms::Exception("Unknown", "FiberSD") << "Cannot find HcalDDDSimConstant\n";
49  }
50 }
51 
53  delete theShower;
54  delete theHC;
55 }
56 
57 void FiberSD::Initialize(G4HCofThisEvent* HCE) {
58  edm::LogVerbatim("FiberSim") << "FiberSD : Initialize called for " << GetName() << " in collection " << HCE;
59  theHC = new FiberG4HitsCollection(GetName(), collectionName[0]);
60  if (theHCID < 0)
61  theHCID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
62  HCE->AddHitsCollection(theHCID, theHC);
63  edm::LogVerbatim("FiberSim") << "FiberSD : Add hit collectrion for " << collectionName[0] << ":" << theHCID << ":"
64  << 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  //edm::LogVerbatim("FiberSim") << "FiberSD :hit position z " << hits[i].position.z();
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  //edm::LogVerbatim("FiberSim") << "FiberSD :presteppoint position z " << preStepPoint->GetPosition().z();
89 
90  FiberG4Hit* aHit = new FiberG4Hit(lv, detID, depth, trackID);
91  edm::LogVerbatim("FiberSim") << "FiberSD :hit size " << hits.size() << " npe" << aHit->npe();
92  edm::LogVerbatim("FiberSim") << "FiberSD :pre hit position " << aHit->hitPos();
93  aHit->setNpe(hits.size());
94  aHit->setPos(theHitPos);
95  aHit->setTime(preStepPoint->GetGlobalTime());
96  aHit->setPhoton(thePE);
97  edm::LogVerbatim("FiberSim") << "FiberSD :ShowerPhoton position " << thePE[0].x() << " " << thePE[0].y() << " "
98  << thePE[0].z();
99 
100  edm::LogVerbatim("FiberSim") << "FiberSD: Hit created at " << lv->GetName() << " DetID: " << aHit->towerId()
101  << " Depth: " << aHit->depth() << " Track ID: " << aHit->trackId()
102  << " Nb. of Cerenkov Photons: " << aHit->npe() << " Time: " << aHit->time() << " at "
103  << aHit->hitPos();
104  for (unsigned int i = 0; i < thePE.size(); i++)
105  edm::LogVerbatim("FiberSim") << "FiberSD: PE[" << i << "] " << thePE[i];
106 
107  theHC->insert(aHit);
108  }
109  return true;
110 }
111 
112 void FiberSD::EndOfEvent(G4HCofThisEvent* HCE) {
113  edm::LogVerbatim("FiberSim") << "FiberSD: finds " << theHC->entries() << " hits";
114  clear();
115  edm::LogVerbatim("FiberSim") << "theHC entries = " << theHC->entries();
116 }
117 
118 void FiberSD::clear() {}
119 
121 
123 
124 void FiberSD::update(const BeginOfJob* job) {}
125 
127 
129 
130 void FiberSD::update(const ::EndOfEvent*) {}
131 
133 
134 uint32_t FiberSD::setDetUnitId(const G4Step* aStep) {
135  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
136  int fibre = (touch->GetReplicaNumber(1)) % 10;
137  int cell = (touch->GetReplicaNumber(2));
138  int tower = (touch->GetReplicaNumber(3));
139  return ((tower * 1000 + cell) * 10 + fibre);
140 }
141 
edm::ESHandle::product
T const * product() const
Definition: ESHandle.h:86
SimTrackManager
Definition: SimTrackManager.h:35
mps_fire.i
i
Definition: mps_fire.py:428
FiberSD::FiberSD
FiberSD(const std::string &, const edm::EventSetup &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
Definition: FiberSD.cc:22
hfClusterShapes_cfi.hits
hits
Definition: hfClusterShapes_cfi.py:5
ESHandle.h
HcalSimNumberingRecord.h
FiberSD::update
void update(const BeginOfJob *) override
This routine will be called when the appropriate signal arrives.
Definition: FiberSD.cc:124
AlCaHLTBitMon_ParallelJobs.p
p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
FiberSD::theHC
FiberG4HitsCollection * theHC
Definition: FiberSD.h:61
FiberSD::~FiberSD
~FiberSD() override
Definition: FiberSD.cc:52
hgcalTowerProducer_cfi.tower
tower
Definition: hgcalTowerProducer_cfi.py:4
FiberSD::EndOfEvent
void EndOfEvent(G4HCofThisEvent *HCE) override
Definition: FiberSD.cc:112
FiberSD::fillHits
void fillHits(edm::PCaloHitContainer &, const std::string &) override
Definition: FiberSD.cc:142
FiberSD::clear
void clear() override
Definition: FiberSD.cc:118
FiberG4Hit
Definition: FiberG4Hit.h:14
HcalDDDSimConstants
Definition: HcalDDDSimConstants.h:24
HcalSimulationConstants::hcalsimpar
const HcalSimulationParameters * hcalsimpar() const
Definition: HcalSimulationConstants.h:20
SensitiveCaloDetector
Definition: SensitiveCaloDetector.h:10
Calorimetry_cff.HFShower
HFShower
Definition: Calorimetry_cff.py:265
edm::EventSetup::get
T get() const
Definition: EventSetup.h:87
FiberSD::DrawAll
void DrawAll() override
Definition: FiberSD.cc:120
FiberSD::theHCID
G4int theHCID
Definition: FiberSD.h:60
HcalSimulationConstants
Definition: HcalSimulationConstants.h:15
edm::ESHandle
Definition: DTSurvey.h:22
SensitiveDetectorCatalog
Definition: SensitiveDetectorCatalog.h:10
BeginOfJob
Definition: BeginOfJob.h:8
bysipixelclustmulteventfilter_cfi.collectionName
collectionName
Definition: bysipixelclustmulteventfilter_cfi.py:5
LEDCalibrationChannels.depth
depth
Definition: LEDCalibrationChannels.py:65
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
edm::ParameterSet
Definition: ParameterSet.h:47
math::XYZPoint
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
FiberSD::ProcessHits
G4bool ProcessHits(G4Step *aStep, G4TouchableHistory *ROhist) override
Definition: FiberSD.cc:67
FiberSD::clearHits
void clearHits() override
Definition: FiberSD.cc:132
position
static int position[264][3]
Definition: ReadPGInfo.cc:289
BeginOfEvent
Definition: BeginOfEvent.h:6
BeginOfRun
Definition: BeginOfRun.h:6
edm::EventSetup
Definition: EventSetup.h:58
edm::LogError
Log< level::Error, false > LogError
Definition: MessageLogger.h:123
get
#define get
edm::ESHandleBase::isValid
bool isValid() const
Definition: ESHandle.h:44
FiberSD::theShower
HFShower * theShower
Definition: FiberSD.h:58
FiberSD.h
edm::LogVerbatim
Log< level::Info, true > LogVerbatim
Definition: MessageLogger.h:128
Exception
Definition: hltDiff.cc:245
HFShowerPhoton
Definition: HFShowerPhoton.h:13
Point3D.h
HcalSimulationConstants.h
edm::PCaloHitContainer
std::vector< PCaloHit > PCaloHitContainer
Definition: PCaloHitContainer.h:8
EventSetup.h
FiberSD::PrintAll
void PrintAll() override
Definition: FiberSD.cc:122
HFShowerPhoton.h
HFShower::getHits
std::vector< Hit > getHits(const G4Step *aStep, double weight)
Definition: HFShower.cc:43
edm::Log
Definition: MessageLogger.h:70
HcalSimNumberingRecord
Definition: HcalSimNumberingRecord.h:25
FiberG4HitsCollection
G4THitsCollection< FiberG4Hit > FiberG4HitsCollection
Definition: FiberG4Hit.h:54
FiberSD::setDetUnitId
uint32_t setDetUnitId(const G4Step *) override
Definition: FiberSD.cc:134
FiberSD::Initialize
void Initialize(G4HCofThisEvent *HCE) override
Definition: FiberSD.cc:57
HcalDDDSimConstants.h