CMS 3D CMS Logo

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