CMS 3D CMS Logo

List of all members | Public Member Functions | Protected Member Functions | Private Attributes
FiberSD Class Reference

#include <FiberSD.h>

Inheritance diagram for FiberSD:
SensitiveCaloDetector Observer< const BeginOfJob *> Observer< const BeginOfRun *> Observer< const BeginOfEvent *> Observer< const EndOfEvent *> SensitiveDetector

Public Member Functions

void clear () override
 
void clearHits () override
 
void DrawAll () override
 
void EndOfEvent (G4HCofThisEvent *HCE) override
 
 FiberSD (const std::string &, const HcalSimulationConstants *, const HcalDDDSimConstants *, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
 
void fillHits (edm::PCaloHitContainer &, const std::string &) override
 
void Initialize (G4HCofThisEvent *HCE) override
 
void PrintAll () override
 
G4bool ProcessHits (G4Step *aStep, G4TouchableHistory *ROhist) override
 
uint32_t setDetUnitId (const G4Step *) override
 
 ~FiberSD () override
 
- Public Member Functions inherited from SensitiveCaloDetector
virtual void reset ()
 
 SensitiveCaloDetector (const std::string &iname, const SensitiveDetectorCatalog &clg, const std::string &newcollname="")
 
- Public Member Functions inherited from SensitiveDetector
void EndOfEvent (G4HCofThisEvent *eventHC) override
 
const std::vector< std::string > & getNames () const
 
void Initialize (G4HCofThisEvent *eventHC) override
 
bool isCaloSD () const
 
 SensitiveDetector (const std::string &iname, const SensitiveDetectorCatalog &, bool calo, const std::string &newcollname="")
 
 ~SensitiveDetector () override
 
- Public Member Functions inherited from Observer< const BeginOfJob *>
 Observer ()
 
void slotForUpdate (const BeginOfJob * iT)
 
virtual ~Observer ()
 
- Public Member Functions inherited from Observer< const BeginOfRun *>
 Observer ()
 
void slotForUpdate (const BeginOfRun * iT)
 
virtual ~Observer ()
 
- Public Member Functions inherited from Observer< const BeginOfEvent *>
 Observer ()
 
void slotForUpdate (const BeginOfEvent * iT)
 
virtual ~Observer ()
 
- Public Member Functions inherited from Observer< const EndOfEvent *>
 Observer ()
 
void slotForUpdate (const EndOfEvent * iT)
 
virtual ~Observer ()
 

Protected Member Functions

void update (const BeginOfJob *) override
 This routine will be called when the appropriate signal arrives. More...
 
void update (const BeginOfRun *) override
 This routine will be called when the appropriate signal arrives. More...
 
void update (const BeginOfEvent *) override
 This routine will be called when the appropriate signal arrives. More...
 
void update (const ::EndOfEvent *) override
 
- Protected Member Functions inherited from SensitiveDetector
TrackInformationcmsTrackInformation (const G4Track *aTrack)
 
Local3DPoint ConvertToLocal3DPoint (const G4ThreeVector &point) const
 
Local3DPoint FinalStepPosition (const G4Step *step, coordinates) const
 
Local3DPoint InitialStepPosition (const G4Step *step, coordinates) const
 
Local3DPoint LocalPostStepPosition (const G4Step *step) const
 
Local3DPoint LocalPreStepPosition (const G4Step *step) const
 
void NaNTrap (const G4Step *step) const
 
void setNames (const std::vector< std::string > &)
 
- Protected Member Functions inherited from Observer< const EndOfEvent *>
virtual void update (const EndOfEvent *)=0
 This routine will be called when the appropriate signal arrives. More...
 

Private Attributes

FiberG4HitsCollectiontheHC
 
G4int theHCID
 
HFShowertheShower
 

Additional Inherited Members

- Protected Types inherited from SensitiveDetector
enum  coordinates { WorldCoordinates, LocalCoordinates }
 

Detailed Description

Definition at line 29 of file FiberSD.h.

Constructor & Destructor Documentation

◆ FiberSD()

FiberSD::FiberSD ( const std::string &  iname,
const HcalSimulationConstants hsps,
const HcalDDDSimConstants hdc,
const SensitiveDetectorCatalog clg,
edm::ParameterSet const &  p,
const SimTrackManager manager 
)
explicit

Definition at line 22 of file FiberSD.cc.

References HcalSimulationConstants::hcalsimpar(), Calorimetry_cff::HFShower, AlCaHLTBitMon_ParallelJobs::p, and theShower.

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 }
Log< level::Info, true > LogVerbatim
SensitiveCaloDetector(const std::string &iname, const SensitiveDetectorCatalog &clg, const std::string &newcollname="")
const HcalSimulationParameters * hcalsimpar() const
HFShower * theShower
Definition: FiberSD.h:61
FiberG4HitsCollection * theHC
Definition: FiberSD.h:64
G4int theHCID
Definition: FiberSD.h:63

◆ ~FiberSD()

FiberSD::~FiberSD ( )
override

Definition at line 34 of file FiberSD.cc.

References theHC, and theShower.

34  {
35  delete theShower;
36  delete theHC;
37 }
HFShower * theShower
Definition: FiberSD.h:61
FiberG4HitsCollection * theHC
Definition: FiberSD.h:64

Member Function Documentation

◆ clear()

void FiberSD::clear ( void  )
override

Definition at line 103 of file FiberSD.cc.

Referenced by EndOfEvent().

103 {}

◆ clearHits()

void FiberSD::clearHits ( )
overridevirtual

Implements SensitiveDetector.

Definition at line 117 of file FiberSD.cc.

117 {}

◆ DrawAll()

void FiberSD::DrawAll ( )
override

Definition at line 105 of file FiberSD.cc.

105 {}

◆ EndOfEvent()

void FiberSD::EndOfEvent ( G4HCofThisEvent *  HCE)
override

Definition at line 97 of file FiberSD.cc.

References clear(), and theHC.

97  {
98  edm::LogVerbatim("FiberSim") << "FiberSD: finds " << theHC->entries() << " hits";
99  clear();
100  edm::LogVerbatim("FiberSim") << "theHC entries = " << theHC->entries();
101 }
Log< level::Info, true > LogVerbatim
FiberG4HitsCollection * theHC
Definition: FiberSD.h:64
void clear() override
Definition: FiberSD.cc:103

◆ fillHits()

void FiberSD::fillHits ( edm::PCaloHitContainer ,
const std::string &   
)
overridevirtual

Implements SensitiveCaloDetector.

Definition at line 127 of file FiberSD.cc.

127 {}

◆ Initialize()

void FiberSD::Initialize ( G4HCofThisEvent *  HCE)
override

Definition at line 39 of file FiberSD.cc.

References bysipixelclustmulteventfilter_cfi::collectionName, theHC, and theHCID.

39  {
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 }
Log< level::Info, true > LogVerbatim
FiberG4HitsCollection * theHC
Definition: FiberSD.h:64
G4THitsCollection< FiberG4Hit > FiberG4HitsCollection
Definition: FiberG4Hit.h:54
G4int theHCID
Definition: FiberSD.h:63

◆ PrintAll()

void FiberSD::PrintAll ( )
override

Definition at line 107 of file FiberSD.cc.

107 {}

◆ ProcessHits()

G4bool FiberSD::ProcessHits ( G4Step *  aStep,
G4TouchableHistory *  ROhist 
)
overridevirtual

Implements SensitiveDetector.

Definition at line 49 of file FiberSD.cc.

References hcalRecHitTable_cff::depth, HFShower::getHits(), hfClusterShapes_cfi::hits, mps_fire::i, position, setDetUnitId(), theHC, and theShower.

49  {
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 }
Log< level::Info, true > LogVerbatim
uint32_t setDetUnitId(const G4Step *) override
Definition: FiberSD.cc:119
std::vector< Hit > getHits(const G4Step *aStep, double weight)
Definition: HFShower.cc:40
HFShower * theShower
Definition: FiberSD.h:61
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
FiberG4HitsCollection * theHC
Definition: FiberSD.h:64
static int position[264][3]
Definition: ReadPGInfo.cc:289

◆ setDetUnitId()

uint32_t FiberSD::setDetUnitId ( const G4Step *  aStep)
overridevirtual

Implements SensitiveDetector.

Definition at line 119 of file FiberSD.cc.

References l1tHGCalTowerProducer_cfi::tower.

Referenced by ProcessHits().

119  {
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 }

◆ update() [1/4]

void FiberSD::update ( const BeginOfJob )
overrideprotectedvirtual

This routine will be called when the appropriate signal arrives.

Implements Observer< const BeginOfJob *>.

Definition at line 109 of file FiberSD.cc.

Referenced by progressbar.ProgressBar::__next__(), MatrixUtil.Matrix::__setitem__(), MatrixUtil.Steps::__setitem__(), progressbar.ProgressBar::finish(), and MatrixUtil.Steps::overwrite().

109 {}

◆ update() [2/4]

void FiberSD::update ( const BeginOfRun )
overrideprotectedvirtual

This routine will be called when the appropriate signal arrives.

Implements Observer< const BeginOfRun *>.

Definition at line 111 of file FiberSD.cc.

Referenced by progressbar.ProgressBar::__next__(), MatrixUtil.Matrix::__setitem__(), MatrixUtil.Steps::__setitem__(), progressbar.ProgressBar::finish(), and MatrixUtil.Steps::overwrite().

111 {}

◆ update() [3/4]

void FiberSD::update ( const BeginOfEvent )
overrideprotectedvirtual

This routine will be called when the appropriate signal arrives.

Implements Observer< const BeginOfEvent *>.

Definition at line 113 of file FiberSD.cc.

Referenced by progressbar.ProgressBar::__next__(), MatrixUtil.Matrix::__setitem__(), MatrixUtil.Steps::__setitem__(), progressbar.ProgressBar::finish(), and MatrixUtil.Steps::overwrite().

113 {}

◆ update() [4/4]

void FiberSD::update ( const ::EndOfEvent )
overrideprotected

Member Data Documentation

◆ theHC

FiberG4HitsCollection* FiberSD::theHC
private

Definition at line 64 of file FiberSD.h.

Referenced by EndOfEvent(), Initialize(), ProcessHits(), and ~FiberSD().

◆ theHCID

G4int FiberSD::theHCID
private

Definition at line 63 of file FiberSD.h.

Referenced by Initialize().

◆ theShower

HFShower* FiberSD::theShower
private

Definition at line 61 of file FiberSD.h.

Referenced by FiberSD(), ProcessHits(), and ~FiberSD().