CMS 3D CMS Logo

TotemRPSD.cc
Go to the documentation of this file.
1 // File: TotemRPSD.cc
2 // Date: 18.10.2005
3 // Description: Sensitive Detector class for TOTEM RP Detectors
4 // Modifications:
7 
10 
15 
17 
20 
21 #include "G4SDManager.hh"
22 #include "G4Step.hh"
23 #include "G4Track.hh"
24 #include "G4VProcess.hh"
25 
26 #include "G4PhysicalConstants.hh"
27 #include "G4SystemOfUnits.hh"
28 
29 #include <iostream>
30 #include <vector>
31 #include <string>
32 
33 static constexpr double rp_garage_position_ = 40.0;
34 
36  const SensitiveDetectorCatalog& clg,
37  edm::ParameterSet const& p,
38  const SimTrackManager* manager)
39  : SensitiveTkDetector(pname, clg) {
40  collectionName.insert(pname);
41 
42  edm::ParameterSet m_Anal = p.getParameter<edm::ParameterSet>("TotemRPSD");
43  verbosity_ = m_Anal.getParameter<int>("Verbosity");
44 
45  slave_ = std::make_unique<TrackingSlaveSD>(pname);
46 
47  if (pname == "TotemHitsRP") {
48  numberingScheme_ = std::make_unique<PPSStripOrganization>();
49  } else {
50  edm::LogWarning("TotemRP") << "TotemRPSD: ReadoutName " << pname << " not supported";
51  }
52 
53  edm::LogVerbatim("TotemRP") << "TotemRPSD: Instantiation completed for " << pname;
54 }
55 
57 
58 void TotemRPSD::Initialize(G4HCofThisEvent* HCE) {
59  LogDebug("TotemRP") << "TotemRPSD : Initialize called for " << GetName();
60 
61  theHC_ = new TotemRPG4HitCollection(GetName(), collectionName[0]);
62  G4SDManager::GetSDMpointer()->AddNewCollection(GetName(), collectionName[0]);
63 
64  if (hcID_ < 0)
65  hcID_ = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
66  HCE->AddHitsCollection(hcID_, theHC_);
67  LogDebug("TotemRP") << "TotemRPSD: is initialized for " << GetName();
68 }
69 
71  LogDebug("TotemRP") << theTrack_->GetDefinition()->GetParticleName() << " TotemRPSD CreateNewHit for"
72  << " PV " << currentPV_->GetName() << " PVid = " << currentPV_->GetCopyNo() << " Unit "
73  << unitID_;
74  LogDebug("TotemRP") << " primary " << primaryID_ << " time slice " << tSliceID_ << " of energy "
75  << theTrack_->GetTotalEnergy() << " Eloss_ " << eloss_ << " position pre: " << hitPoint_
76  << " post: " << exitPoint_;
77 
78  if (parentID_ == 0) {
79  LogDebug("TotemRP") << " primary particle ";
80  } else {
81  LogDebug("TotemRP") << " daughter of part. " << parentID_;
82  }
83 
84  LogDebug("TotemRP") << " and created by ";
85 
86  if (theTrack_->GetCreatorProcess() != nullptr)
87  LogDebug("TotemRP") << theTrack_->GetCreatorProcess()->GetProcessName();
88  else
89  LogDebug("TotemRP") << "NO process";
90 }
91 
92 bool TotemRPSD::ProcessHits(G4Step* aStep, G4TouchableHistory*) {
93  eloss_ = aStep->GetTotalEnergyDeposit();
94  if (eloss_ > 0.0) {
95  eloss_ /= GeV;
96  stepInfo(aStep);
97  edm::LogVerbatim("TotemRP") << "TotemRPSD: ProcessHits 1: Eloss=" << eloss_ << " "
98  << theTrack_->GetDefinition()->GetParticleName();
99  createNewHit();
100  }
101  return true;
102 }
103 
104 void TotemRPSD::stepInfo(const G4Step* aStep) {
105  preStepPoint_ = aStep->GetPreStepPoint();
106  postStepPoint_ = aStep->GetPostStepPoint();
107  theTrack_ = aStep->GetTrack();
108  hitPoint_ = preStepPoint_->GetPosition();
109  exitPoint_ = postStepPoint_->GetPosition();
110  currentPV_ = preStepPoint_->GetPhysicalVolume();
113 
114  tSlice_ = (postStepPoint_->GetGlobalTime()) / nanosecond;
115  tSliceID_ = (int)tSlice_;
116  unitID_ = setDetUnitId(aStep);
117 
118  if (verbosity_)
119  LogDebug("TotemRP") << "UNIT " << unitID_;
120 
121  primaryID_ = theTrack_->GetTrackID();
122  parentID_ = theTrack_->GetParentID();
123 
124  incidentEnergy_ = theTrack_->GetTotalEnergy() / GeV;
125 
126  pabs_ = preStepPoint_->GetMomentum().mag() / GeV;
127  thePx_ = preStepPoint_->GetMomentum().x() / GeV;
128  thePy_ = preStepPoint_->GetMomentum().y() / GeV;
129  thePz_ = preStepPoint_->GetMomentum().z() / GeV;
130 
131  tof_ = postStepPoint_->GetGlobalTime() / nanosecond;
132  particleType_ = theTrack_->GetDefinition()->GetPDGEncoding();
133 
134  //corrected phi and theta treatment
135  G4ThreeVector gmd = preStepPoint_->GetMomentumDirection();
136  // convert it to local frame
137  G4ThreeVector lmd =
138  ((G4TouchableHistory*)(preStepPoint_->GetTouchable()))->GetHistory()->GetTopTransform().TransformAxis(gmd);
140  thetaAtEntry_ = lnmd.theta();
141  phiAtEntry_ = lnmd.phi();
142 
143  vx_ = theTrack_->GetVertexPosition().x() / mm;
144  vy_ = theTrack_->GetVertexPosition().y() / mm;
145  vz_ = theTrack_->GetVertexPosition().z() / mm;
146 }
147 
148 uint32_t TotemRPSD::setDetUnitId(const G4Step* aStep) {
149  return (numberingScheme_ == nullptr ? 0 : numberingScheme_->unitID(aStep));
150 }
151 
153  if (hit == nullptr) {
154  if (verbosity_)
155  LogDebug("TotemRP") << "TotemRPSD: hit to be stored is NULL !!" << std::endl;
156  } else {
157  theHC_->insert(hit);
158  }
159 }
160 
162  // Protect against creating hits in detectors not inserted
163  double outrangeX = hitPoint_.x();
164  double outrangeY = hitPoint_.y();
165  if (std::abs(outrangeX) > rp_garage_position_ || std::abs(outrangeY) > rp_garage_position_)
166  return;
167  // end protection
168 
174 
181 
186 
191 
195 
197 }
198 
199 G4ThreeVector TotemRPSD::setToLocal(const G4ThreeVector& global) {
200  G4ThreeVector localPoint;
201  const G4VTouchable* touch = preStepPoint_->GetTouchable();
202  localPoint = touch->GetHistory()->GetTopTransform().TransformPoint(global);
203  return localPoint;
204 }
205 
206 void TotemRPSD::EndOfEvent(G4HCofThisEvent*) {
207  // here we loop over transient hits and make them persistent
208  for (unsigned int j = 0; j < (unsigned int)theHC_->entries(); ++j) {
209  TotemRPG4Hit* aHit = (*theHC_)[j];
210 
211  Local3DPoint entry(aHit->localEntry().x(), aHit->localEntry().y(), aHit->localEntry().z());
212  Local3DPoint exit(aHit->localExit().x(), aHit->localExit().y(), aHit->localExit().z());
213  slave_->processHits(PSimHit(entry,
214  exit,
215  aHit->p(),
216  aHit->tof(),
217  aHit->energyLoss(),
218  aHit->particleType(),
219  aHit->unitID(),
220  aHit->trackID(),
221  aHit->thetaAtEntry(),
222  aHit->phiAtEntry()));
223  }
224 }
225 
227  LogDebug("TotemRP") << "TotemRPSD: Collection " << theHC_->GetName() << std::endl;
228  theHC_->PrintAllHits();
229 }
230 
232  if (slave_->name() == n) {
233  c = slave_->hits();
234  }
235 }
236 
237 void TotemRPSD::update(const BeginOfEvent* ptr) {
238  clearHits();
239  eventno_ = (*ptr)()->GetEventID();
240 }
241 
242 void TotemRPSD::update(const ::EndOfEvent*) {}
243 
244 void TotemRPSD::clearHits() { slave_->Initialize(); }
void setEntry(G4ThreeVector xyz)
Definition: TotemRPG4Hit.cc:76
G4ThreeVector localExit() const
Definition: TotemRPG4Hit.cc:83
Log< level::Info, true > LogVerbatim
void fillHits(edm::PSimHitContainer &, const std::string &) override
Definition: TotemRPSD.cc:231
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
G4ThreeVector localEntry() const
Definition: TotemRPG4Hit.cc:81
G4Track * theTrack_
Definition: TotemRPSD.h:63
int unitID() const
Definition: TotemRPG4Hit.cc:92
void setTof(double e)
G4THitsCollection< TotemRPG4Hit > TotemRPG4HitCollection
double pabs_
Definition: TotemRPSD.h:79
void setP(double e)
G4ThreeVector setToLocal(const G4ThreeVector &globalPoint)
Definition: TotemRPSD.cc:199
void createNewHit()
Definition: TotemRPSD.cc:161
void setIncidentEnergy(double e)
Definition: TotemRPG4Hit.cc:87
G4int tSliceID_
Definition: TotemRPSD.h:68
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
G4VPhysicalVolume * currentPV_
Definition: TotemRPSD.h:64
std::unique_ptr< TotemRPVDetectorOrganization > numberingScheme_
Definition: TotemRPSD.h:59
void storeHit(TotemRPG4Hit *)
Definition: TotemRPSD.cc:152
bool ProcessHits(G4Step *step, G4TouchableHistory *tHistory) override
Definition: TotemRPSD.cc:92
double vx_
Definition: TotemRPSD.h:89
void PrintAll() override
Definition: TotemRPSD.cc:226
int eventno_
Definition: TotemRPSD.h:95
void setVx(double p)
G4int parentID_
Definition: TotemRPSD.h:67
double energyLoss() const
void setPy(double p)
TotemRPG4HitCollection * theHC_
Definition: TotemRPSD.h:61
int particleType() const
double thetaAtEntry() const
G4double tSlice_
Definition: TotemRPSD.h:69
double phiAtEntry() const
void setLocalEntry(const G4ThreeVector &theLocalEntryPoint)
Definition: TotemRPG4Hit.cc:82
void stepInfo(const G4Step *aStep)
Definition: TotemRPSD.cc:104
void clearHits() override
Definition: TotemRPSD.cc:244
void printHitInfo()
Definition: TotemRPSD.cc:70
void Initialize(G4HCofThisEvent *HCE) override
Definition: TotemRPSD.cc:58
void EndOfEvent(G4HCofThisEvent *eventHC) override
Definition: TotemRPSD.cc:206
void setExit(G4ThreeVector xyz)
Definition: TotemRPG4Hit.cc:79
double vz_
Definition: TotemRPSD.h:91
Local3DPoint ConvertToLocal3DPoint(const G4ThreeVector &point) const
G4ThreeVector theLocalEntryPoint_
Definition: TotemRPSD.h:75
uint32_t setDetUnitId(const G4Step *step) override
Definition: TotemRPSD.cc:148
G4StepPoint * postStepPoint_
Definition: TotemRPSD.h:72
int verbosity_
Definition: TotemRPSD.h:94
double vy_
Definition: TotemRPSD.h:90
unsigned int unitID_
Definition: TotemRPSD.h:93
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
G4int primaryID_
Definition: TotemRPSD.h:66
double tof() const
void setLocalExit(const G4ThreeVector &theLocalExitPoint)
Definition: TotemRPG4Hit.cc:84
double thetaAtEntry_
Definition: TotemRPSD.h:86
void setThetaAtEntry(double t)
G4StepPoint * preStepPoint_
Definition: TotemRPSD.h:71
void setUnitID(unsigned int i)
Definition: TotemRPG4Hit.cc:93
G4ThreeVector theLocalExitPoint_
Definition: TotemRPSD.h:76
void setPhiAtEntry(double f)
double p() const
Definition: TotemRPG4Hit.cc:99
double thePz_
Definition: TotemRPSD.h:82
void setVz(double p)
double thePx_
Definition: TotemRPSD.h:80
double incidentEnergy_
Definition: TotemRPSD.h:78
double eloss_
Definition: TotemRPSD.h:84
void setParentId(int p)
TotemRPSD(const std::string &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
Definition: TotemRPSD.cc:35
void setTrackID(int i)
Definition: TotemRPG4Hit.cc:90
G4ThreeVector exitPoint_
Definition: TotemRPSD.h:74
double tof_
Definition: TotemRPSD.h:83
void setPz(double p)
void update(const BeginOfEvent *) override
This routine will be called when the appropriate signal arrives.
Definition: TotemRPSD.cc:237
static constexpr double rp_garage_position_
Definition: TotemRPSD.cc:33
void setParticleType(short i)
G4int hcID_
Definition: TotemRPSD.h:65
std::vector< PSimHit > PSimHitContainer
short particleType_
Definition: TotemRPSD.h:96
Log< level::Warning, false > LogWarning
~TotemRPSD() override
Definition: TotemRPSD.cc:56
void setEnergyLoss(double e)
double phiAtEntry_
Definition: TotemRPSD.h:87
std::unique_ptr< TrackingSlaveSD > slave_
Definition: TotemRPSD.h:58
void setVy(double p)
double thePy_
Definition: TotemRPSD.h:81
void setPx(double p)
G4ThreeVector hitPoint_
Definition: TotemRPSD.h:73
TotemRPG4Hit * currentHit_
Definition: TotemRPSD.h:62
Geom::Theta< T > theta() const
Definition: PV3DBase.h:72
unsigned int trackID() const
Definition: TotemRPG4Hit.cc:89
void setTimeSlice(double d)
Definition: TotemRPG4Hit.cc:96
#define LogDebug(id)
def exit(msg="")