CMS 3D CMS Logo

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

#include <HFChamberSD.h>

Inheritance diagram for HFChamberSD:
SensitiveCaloDetector SensitiveDetector

Public Member Functions

void clear () override
 
void clearHits () override
 
void DrawAll () override
 
void EndOfEvent (G4HCofThisEvent *HCE) override
 
void fillHits (edm::PCaloHitContainer &, const std::string &) override
 
 HFChamberSD (const std::string &, const DDCompactView &, const SensitiveDetectorCatalog &, const edm::ParameterSet &, const SimTrackManager *)
 
void Initialize (G4HCofThisEvent *HCE) override
 
void PrintAll () override
 
G4bool ProcessHits (G4Step *aStep, G4TouchableHistory *ROhist) override
 
uint32_t setDetUnitId (const G4Step *) override
 
 ~HFChamberSD () override
 
- Public Member Functions inherited from SensitiveCaloDetector
 SensitiveCaloDetector (const std::string &iname, const DDCompactView &cpv, const SensitiveDetectorCatalog &clg, edm::ParameterSet const &p)
 
- Public Member Functions inherited from SensitiveDetector
void EndOfEvent (G4HCofThisEvent *eventHC) override
 
const std::vector< std::string > & getNames () const
 
void Initialize (G4HCofThisEvent *eventHC) override
 
 SensitiveDetector (const std::string &iname, const DDCompactView &cpv, const SensitiveDetectorCatalog &, edm::ParameterSet const &p)
 
 ~SensitiveDetector () override
 

Private Attributes

const SimTrackManagerm_trackManager
 
HFShowerG4HitsCollectiontheHC
 
G4int theHCID
 
int theNSteps
 

Additional Inherited Members

- Protected Types inherited from SensitiveDetector
enum  coordinates { WorldCoordinates, LocalCoordinates }
 
- Protected Member Functions inherited from SensitiveDetector
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 > &)
 

Detailed Description

Definition at line 19 of file HFChamberSD.h.

Constructor & Destructor Documentation

HFChamberSD::HFChamberSD ( const std::string &  ,
const DDCompactView ,
const SensitiveDetectorCatalog ,
const edm::ParameterSet ,
const SimTrackManager  
)
explicit

Definition at line 19 of file HFChamberSD.cc.

21  :
22  SensitiveCaloDetector(name, cpv, clg, p),
23  m_trackManager(manager), theHCID(-1), theHC(nullptr), theNSteps(0) {
24 
25 }
SensitiveCaloDetector(const std::string &iname, const DDCompactView &cpv, const SensitiveDetectorCatalog &clg, edm::ParameterSet const &p)
HFShowerG4HitsCollection * theHC
Definition: HFChamberSD.h:44
G4int theHCID
Definition: HFChamberSD.h:43
const SimTrackManager * m_trackManager
Definition: HFChamberSD.h:41
HFChamberSD::~HFChamberSD ( )
override

Definition at line 27 of file HFChamberSD.cc.

References theHC.

27  {
28  delete theHC;
29 }
HFShowerG4HitsCollection * theHC
Definition: HFChamberSD.h:44

Member Function Documentation

void HFChamberSD::clear ( void  )
override

Definition at line 84 of file HFChamberSD.cc.

References theNSteps.

Referenced by EndOfEvent().

84  {
85  theNSteps = 0;
86 }
void HFChamberSD::clearHits ( )
overridevirtual

Implements SensitiveDetector.

Definition at line 92 of file HFChamberSD.cc.

92 {}
void HFChamberSD::DrawAll ( )
override

Definition at line 88 of file HFChamberSD.cc.

88 {}
void HFChamberSD::EndOfEvent ( G4HCofThisEvent *  HCE)
override

Definition at line 78 of file HFChamberSD.cc.

References clear(), LogDebug, and theHC.

78  {
79 
80  LogDebug("FiberSim") << "HFChamberSD: Sees" << theHC->entries() << " hits";
81  clear();
82 }
#define LogDebug(id)
void clear() override
Definition: HFChamberSD.cc:84
HFShowerG4HitsCollection * theHC
Definition: HFChamberSD.h:44
void HFChamberSD::fillHits ( edm::PCaloHitContainer ,
const std::string &   
)
overridevirtual

Implements SensitiveCaloDetector.

Definition at line 99 of file HFChamberSD.cc.

99 {}
void HFChamberSD::Initialize ( G4HCofThisEvent *  HCE)
override

Definition at line 31 of file HFChamberSD.cc.

References ecaldqm::collectionName, LogDebug, theHC, and theHCID.

31  {
32 
33  LogDebug("FiberSim") << "HFChamberSD : Initialize called for " << GetName();
34  theHC = new HFShowerG4HitsCollection(GetName(), collectionName[0]);
35  if (theHCID<0)
36  theHCID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
37  HCE->AddHitsCollection(theHCID, theHC);
38 
39 }
#define LogDebug(id)
HFShowerG4HitsCollection * theHC
Definition: HFChamberSD.h:44
std::string const collectionName[nCollections]
Definition: Collections.h:47
G4int theHCID
Definition: HFChamberSD.h:43
G4THitsCollection< HFShowerG4Hit > HFShowerG4HitsCollection
Definition: HFShowerG4Hit.h:59
void HFChamberSD::PrintAll ( )
override

Definition at line 90 of file HFChamberSD.cc.

90 {}
G4bool HFChamberSD::ProcessHits ( G4Step *  aStep,
G4TouchableHistory *  ROhist 
)
overridevirtual

Implements SensitiveDetector.

Definition at line 41 of file HFChamberSD.cc.

References ALCARECOTkAlJpsiMuMu_cff::charge, LogDebug, MeV, setDetUnitId(), HFShowerG4Hit::setGlobalPos(), HFShowerG4Hit::setLocalPos(), HFShowerG4Hit::setPrimMomDir(), theHC, theNSteps, and ntuplemaker::time.

41  {
42 
43  //do not process hits other than primary particle hits:
44  double charge = aStep->GetTrack()->GetDefinition()->GetPDGCharge();
45  int trackID = aStep->GetTrack()->GetTrackID();
46  if(charge == 0. || trackID != 1 ||aStep->GetTrack()->GetParentID() != 0 || aStep->GetTrack()->GetCreatorProcess() != nullptr) return false;
47  ++theNSteps;
48  //if(theNSteps>1)return false;
49 
50  G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
51  const G4VTouchable* touch = preStepPoint->GetTouchable();
52  int detID = setDetUnitId(aStep);
53 
54  double edep = aStep->GetTotalEnergyDeposit();
55  double time = (preStepPoint->GetGlobalTime())/ns;
56 
57  const G4ThreeVector& globalPos = preStepPoint->GetPosition();
58  G4ThreeVector localPos = touch->GetHistory()->GetTopTransform().TransformPoint(globalPos);
59  const G4DynamicParticle* particle = aStep->GetTrack()->GetDynamicParticle();
60  const G4ThreeVector& momDir = particle->GetMomentumDirection();
61 
62  HFShowerG4Hit *aHit = new HFShowerG4Hit(detID, trackID, edep, time);
63  aHit->setLocalPos(localPos);
64  aHit->setGlobalPos(globalPos);
65  aHit->setPrimMomDir(momDir);
66 
67  LogDebug("FiberSim") << "HFChamberSD: Hit created in ("
68  << touch->GetVolume(0)->GetLogicalVolume()->GetName()
69  << ") " << " ID " << detID << " Track " << trackID
70  << " Edep: " << edep/MeV << " MeV; Time: " << time
71  << " ns; Position (local) " << localPos << " (global ) "
72  << globalPos << " direction " << momDir;
73 
74  theHC->insert(aHit);
75  return true;
76 }
#define LogDebug(id)
void setLocalPos(const G4ThreeVector &xyz)
Definition: HFShowerG4Hit.h:46
uint32_t setDetUnitId(const G4Step *) override
Definition: HFChamberSD.cc:94
const double MeV
HFShowerG4HitsCollection * theHC
Definition: HFChamberSD.h:44
void setPrimMomDir(const G4ThreeVector &xyz)
Definition: HFShowerG4Hit.h:48
void setGlobalPos(const G4ThreeVector &xyz)
Definition: HFShowerG4Hit.h:47
uint32_t HFChamberSD::setDetUnitId ( const G4Step *  aStep)
overridevirtual

Implements SensitiveDetector.

Definition at line 94 of file HFChamberSD.cc.

Referenced by ProcessHits().

94  {
95  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
96  return (touch->GetReplicaNumber(0));
97 }

Member Data Documentation

const SimTrackManager* HFChamberSD::m_trackManager
private

Definition at line 41 of file HFChamberSD.h.

HFShowerG4HitsCollection* HFChamberSD::theHC
private

Definition at line 44 of file HFChamberSD.h.

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

G4int HFChamberSD::theHCID
private

Definition at line 43 of file HFChamberSD.h.

Referenced by Initialize().

int HFChamberSD::theNSteps
private

Definition at line 45 of file HFChamberSD.h.

Referenced by clear(), and ProcessHits().