#include <HFChamberSD.h>
Public Member Functions | |
virtual void | clear () |
virtual void | DrawAll () |
virtual void | EndOfEvent (G4HCofThisEvent *HCE) |
HFChamberSD (std::string, const DDCompactView &, SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *) | |
virtual void | Initialize (G4HCofThisEvent *HCE) |
virtual void | PrintAll () |
virtual G4bool | ProcessHits (G4Step *aStep, G4TouchableHistory *ROhist) |
virtual | ~HFChamberSD () |
Protected Member Functions | |
virtual void | clearHits () |
virtual void | fillHits (edm::PCaloHitContainer &, std::string) |
virtual uint32_t | setDetUnitId (G4Step *) |
Private Attributes | |
const SimTrackManager * | m_trackManager |
HFShowerG4HitsCollection * | theHC |
G4int | theHCID |
std::string | theName |
int | theNSteps |
Definition at line 19 of file HFChamberSD.h.
HFChamberSD::HFChamberSD | ( | std::string | name, |
const DDCompactView & | cpv, | ||
SensitiveDetectorCatalog & | clg, | ||
edm::ParameterSet const & | p, | ||
const SimTrackManager * | manager | ||
) |
Definition at line 16 of file HFChamberSD.cc.
References SensitiveDetector::AssignSD(), LogDebug, SensitiveDetectorCatalog::logicalNames(), and SensitiveDetector::Register().
: SensitiveCaloDetector(name, cpv, clg, p), theName(name), m_trackManager(manager), theHCID(-1), theHC(0), theNSteps(0) { collectionName.insert(name); LogDebug("FiberSim") << "***************************************************" << "\n" << "* *" << "\n" << "* Constructing a HFChamberSD with name " << GetName() << "\n" << "* *" << "\n" << "***************************************************"; // // Now attach the right detectors (LogicalVolumes) to me // std::vector<std::string> lvNames = clg.logicalNames(name); this->Register(); for (std::vector<std::string>::iterator it=lvNames.begin(); it !=lvNames.end(); it++){ this->AssignSD(*it); LogDebug("FiberSim") << "HFChamberSD : Assigns SD to LV " << (*it); } }
HFChamberSD::~HFChamberSD | ( | ) | [virtual] |
void HFChamberSD::clear | ( | void | ) | [virtual] |
Definition at line 101 of file HFChamberSD.cc.
References theNSteps.
Referenced by EndOfEvent().
{ theNSteps = 0; }
void HFChamberSD::clearHits | ( | ) | [protected, virtual] |
void HFChamberSD::DrawAll | ( | ) | [virtual] |
Definition at line 105 of file HFChamberSD.cc.
{}
void HFChamberSD::EndOfEvent | ( | G4HCofThisEvent * | HCE | ) | [virtual] |
Reimplemented from SensitiveDetector.
Definition at line 95 of file HFChamberSD.cc.
void HFChamberSD::fillHits | ( | edm::PCaloHitContainer & | , |
std::string | |||
) | [protected, virtual] |
void HFChamberSD::Initialize | ( | G4HCofThisEvent * | HCE | ) | [virtual] |
Reimplemented from SensitiveDetector.
Definition at line 48 of file HFChamberSD.cc.
void HFChamberSD::PrintAll | ( | ) | [virtual] |
Definition at line 107 of file HFChamberSD.cc.
{}
G4bool HFChamberSD::ProcessHits | ( | G4Step * | aStep, |
G4TouchableHistory * | ROhist | ||
) | [virtual] |
Implements SensitiveDetector.
Definition at line 58 of file HFChamberSD.cc.
References DeDxDiscriminatorTools::charge(), LogDebug, NULL, setDetUnitId(), HFShowerG4Hit::setGlobalPos(), HFShowerG4Hit::setLocalPos(), HFShowerG4Hit::setPrimMomDir(), theHC, theNSteps, and cond::rpcobgas::time.
{ //do not process hits other than primary particle hits: double charge = aStep->GetTrack()->GetDefinition()->GetPDGCharge(); int trackID = aStep->GetTrack()->GetTrackID(); if(charge == 0. || trackID != 1 ||aStep->GetTrack()->GetParentID() != 0 || aStep->GetTrack()->GetCreatorProcess() != NULL) return false; ++theNSteps; //if(theNSteps>1)return false; G4StepPoint* preStepPoint = aStep->GetPreStepPoint(); const G4VTouchable* touch = preStepPoint->GetTouchable(); int detID = setDetUnitId(aStep); double edep = aStep->GetTotalEnergyDeposit(); double time = (preStepPoint->GetGlobalTime())/ns; G4ThreeVector globalPos = preStepPoint->GetPosition(); G4ThreeVector localPos = touch->GetHistory()->GetTopTransform().TransformPoint(globalPos); const G4DynamicParticle* particle = aStep->GetTrack()->GetDynamicParticle(); G4ThreeVector momDir = particle->GetMomentumDirection(); HFShowerG4Hit *aHit = new HFShowerG4Hit(detID, trackID, edep, time); aHit->setLocalPos(localPos); aHit->setGlobalPos(globalPos); aHit->setPrimMomDir(momDir); LogDebug("FiberSim") << "HFChamberSD: Hit created in (" << touch->GetVolume(0)->GetLogicalVolume()->GetName() << ") " << " ID " << detID << " Track " << trackID << " Edep: " << edep/MeV << " MeV; Time: " << time << " ns; Position (local) " << localPos << " (global ) " << globalPos << " direction " << momDir; theHC->insert(aHit); return true; }
uint32_t HFChamberSD::setDetUnitId | ( | G4Step * | aStep | ) | [protected, virtual] |
Implements SensitiveDetector.
Definition at line 111 of file HFChamberSD.cc.
Referenced by ProcessHits().
{ const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable(); return (touch->GetReplicaNumber(0)); }
const SimTrackManager* HFChamberSD::m_trackManager [private] |
Definition at line 43 of file HFChamberSD.h.
HFShowerG4HitsCollection* HFChamberSD::theHC [private] |
Definition at line 46 of file HFChamberSD.h.
Referenced by EndOfEvent(), Initialize(), ProcessHits(), and ~HFChamberSD().
G4int HFChamberSD::theHCID [private] |
Definition at line 45 of file HFChamberSD.h.
Referenced by Initialize().
std::string HFChamberSD::theName [private] |
Definition at line 42 of file HFChamberSD.h.
int HFChamberSD::theNSteps [private] |
Definition at line 47 of file HFChamberSD.h.
Referenced by clear(), and ProcessHits().