CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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 
34  const SensitiveDetectorCatalog& clg,
35  edm::ParameterSet const& p,
36  const SimTrackManager* manager)
37  : SensitiveTkDetector(name_, clg),
38  numberingScheme_(nullptr),
39  hcID_(-1),
40  theHC_(nullptr),
41  currentHit_(nullptr),
42  theTrack_(nullptr),
43  currentPV_(nullptr),
44  unitID_(0),
45  preStepPoint_(nullptr),
46  postStepPoint_(nullptr),
47  eventno_(0) {
48  collectionName.insert(name_);
49 
50  edm::ParameterSet m_Anal = p.getParameter<edm::ParameterSet>("TotemRPSD");
51  verbosity_ = m_Anal.getParameter<int>("Verbosity");
52 
53  slave_ = std::make_unique<TrackingSlaveSD>(name_);
54 
55  if (name_ == "TotemHitsRP") {
56  numberingScheme_ = std::make_unique<PPSStripNumberingScheme>(3);
57  } else {
58  edm::LogWarning("TotemRP") << "TotemRPSD: ReadoutName not supported\n";
59  }
60 
61  edm::LogInfo("TotemRP") << "TotemRPSD: Instantiation completed";
62 }
63 
65 
66 void TotemRPSD::Initialize(G4HCofThisEvent* HCE) {
67  LogDebug("TotemRP") << "TotemRPSD : Initialize called for " << name_;
68 
69  theHC_ = new TotemRPG4HitCollection(GetName(), collectionName[0]);
70  G4SDManager::GetSDMpointer()->AddNewCollection(name_, collectionName[0]);
71 
72  if (hcID_ < 0)
73  hcID_ = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
74  HCE->AddHitsCollection(hcID_, theHC_);
75 }
76 
78  LogDebug("TotemRP") << theTrack_->GetDefinition()->GetParticleName() << " TotemRPSD CreateNewHit for"
79  << " PV " << currentPV_->GetName() << " PVid = " << currentPV_->GetCopyNo() << " Unit "
80  << unitID_;
81  LogDebug("TotemRP") << " primary " << primaryID_ << " time slice " << tSliceID_ << " of energy "
82  << theTrack_->GetTotalEnergy() << " Eloss_ " << Eloss_ << " positions ";
83  printf("(%10f,%10f,%10f)",
84  preStepPoint_->GetPosition().x(),
85  preStepPoint_->GetPosition().y(),
86  preStepPoint_->GetPosition().z());
87  printf("(%10f,%10f,%10f)",
88  postStepPoint_->GetPosition().x(),
89  postStepPoint_->GetPosition().y(),
90  postStepPoint_->GetPosition().z());
91  LogDebug("TotemRP") << " positions "
92  << "(" << postStepPoint_->GetPosition().x() << "," << postStepPoint_->GetPosition().y() << ","
93  << postStepPoint_->GetPosition().z() << ")"
94  << " For Track " << theTrack_->GetTrackID() << " which is a "
95  << theTrack_->GetDefinition()->GetParticleName();
96 
97  if (theTrack_->GetTrackID() == 1) {
98  LogDebug("TotemRP") << " primary particle ";
99  } else {
100  LogDebug("TotemRP") << " daughter of part. " << theTrack_->GetParentID();
101  }
102 
103  LogDebug("TotemRP") << " and created by ";
104 
105  if (theTrack_->GetCreatorProcess() != nullptr)
106  LogDebug("TotemRP") << theTrack_->GetCreatorProcess()->GetProcessName();
107  else
108  LogDebug("TotemRP") << "NO process";
109 
110  LogDebug("TotemRP") << std::endl;
111 }
112 
113 bool TotemRPSD::ProcessHits(G4Step* aStep, G4TouchableHistory*) {
114  if (aStep == nullptr) {
115  return true;
116  } else {
117  stepInfo(aStep);
118 
119  createNewHit();
120  return true;
121  }
122 }
123 
124 void TotemRPSD::stepInfo(const G4Step* aStep) {
125  preStepPoint_ = aStep->GetPreStepPoint();
126  postStepPoint_ = aStep->GetPostStepPoint();
127  theTrack_ = aStep->GetTrack();
128  hitPoint_ = preStepPoint_->GetPosition();
129  exitPoint_ = postStepPoint_->GetPosition();
130  currentPV_ = preStepPoint_->GetPhysicalVolume();
133 
134  G4String name_ = currentPV_->GetName();
135  name_.assign(name_, 0, 4);
136  G4String particleType = theTrack_->GetDefinition()->GetParticleName();
137  tSlice_ = (postStepPoint_->GetGlobalTime()) / nanosecond;
138  tSliceID_ = (int)tSlice_;
139  unitID_ = setDetUnitId(aStep);
140 
141  if (verbosity_)
142  LogDebug("TotemRP") << "UNIT " << unitID_ << std::endl;
143 
144  primaryID_ = theTrack_->GetTrackID();
145 
146  Pabs_ = (aStep->GetPreStepPoint()->GetMomentum().mag()) / GeV;
147  thePx_ = (aStep->GetPreStepPoint()->GetMomentum().x()) / GeV;
148  thePy_ = (aStep->GetPreStepPoint()->GetMomentum().y()) / GeV;
149  thePz_ = (aStep->GetPreStepPoint()->GetMomentum().z()) / GeV;
150 
151  Tof_ = aStep->GetPostStepPoint()->GetGlobalTime() / nanosecond;
152  Eloss_ = aStep->GetTotalEnergyDeposit() / GeV;
153  ParticleType_ = theTrack_->GetDefinition()->GetPDGEncoding();
154 
155  //corrected phi and theta treatment
156  G4ThreeVector gmd = aStep->GetPreStepPoint()->GetMomentumDirection();
157  // convert it to local frame
158  G4ThreeVector lmd = ((G4TouchableHistory*)(aStep->GetPreStepPoint()->GetTouchable()))
159  ->GetHistory()
160  ->GetTopTransform()
161  .TransformAxis(gmd);
163  ThetaAtEntry_ = lnmd.theta();
164  PhiAtEntry_ = lnmd.phi();
165 
166  if (isPrimary(theTrack_))
167  ParentId_ = 0;
168  else
169  ParentId_ = theTrack_->GetParentID();
170 
171  Vx_ = theTrack_->GetVertexPosition().x() / mm;
172  Vy_ = theTrack_->GetVertexPosition().y() / mm;
173  Vz_ = theTrack_->GetVertexPosition().z() / mm;
174 }
175 
176 uint32_t TotemRPSD::setDetUnitId(const G4Step* aStep) {
177  return (numberingScheme_ == nullptr ? 0 : numberingScheme_->unitID(aStep));
178 }
179 
181  if (hit == nullptr) {
182  if (verbosity_)
183  LogDebug("TotemRP") << "TotemRPSD: hit to be stored is NULL !!" << std::endl;
184  return;
185  }
186  theHC_->insert(hit);
187 }
188 
190  // Protect against creating hits in detectors not inserted
191  double outrangeX = hitPoint_.x();
192  double outrangeY = hitPoint_.y();
193  if (fabs(outrangeX) > rp_garage_position_)
194  return;
195  if (fabs(outrangeY) > rp_garage_position_)
196  return;
197  // end protection
198 
204 
211 
216 
221 
225 
227 }
228 
229 G4ThreeVector TotemRPSD::setToLocal(const G4ThreeVector& global) {
230  G4ThreeVector localPoint;
231  const G4VTouchable* touch = preStepPoint_->GetTouchable();
232  localPoint = touch->GetHistory()->GetTopTransform().TransformPoint(global);
233  return localPoint;
234 }
235 
236 void TotemRPSD::EndOfEvent(G4HCofThisEvent*) {
237  // here we loop over transient hits and make them persistent
238  for (unsigned int j = 0; j < (unsigned int)theHC_->entries() && j < maxTotemHits_; j++) {
239  TotemRPG4Hit* aHit = (*theHC_)[j];
240 
241  Local3DPoint entry(aHit->localEntry().x(), aHit->localEntry().y(), aHit->localEntry().z());
242  Local3DPoint exit(aHit->localExit().x(), aHit->localExit().y(), aHit->localExit().z());
243  slave_->processHits(PSimHit(entry,
244  exit,
245  aHit->p(),
246  aHit->tof(),
247  aHit->energyLoss(),
248  aHit->particleType(),
249  aHit->unitID(),
250  aHit->trackID(),
251  aHit->thetaAtEntry(),
252  aHit->phiAtEntry()));
253  }
254  summarize();
255 }
256 
258 
260 
262 
264  LogDebug("TotemRP") << "TotemRPSD: Collection " << theHC_->GetName() << std::endl;
265  theHC_->PrintAllHits();
266 }
267 
269  if (slave_->name() == n) {
270  c = slave_->hits();
271  }
272 }
274  if (scheme) {
275  LogDebug("TotemRP") << "TotemRPSD: updates numbering scheme for " << GetName();
276  numberingScheme_.reset(scheme);
277  }
278 }
280  clearHits();
281  eventno_ = (*i)()->GetEventID();
282 }
283 
284 void TotemRPSD::update(const ::EndOfEvent*) {}
285 
286 void TotemRPSD::clearHits() { slave_->Initialize(); }
287 
288 bool TotemRPSD::isPrimary(const G4Track* track) {
289  TrackInformation* info = dynamic_cast<TrackInformation*>(track->GetUserInformation());
290  return info && info->isPrimary();
291 }
void setEntry(G4ThreeVector xyz)
Definition: TotemRPG4Hit.cc:76
G4ThreeVector localEntry() const
Definition: TotemRPG4Hit.cc:81
void fillHits(edm::PSimHitContainer &, const std::string &) override
Definition: TotemRPSD.cc:268
unsigned int trackID() const
Definition: TotemRPG4Hit.cc:89
int particleType() const
static const TGPicture * info(bool iBackgroundIsBlack)
const edm::EventSetup & c
G4Track * theTrack_
Definition: TotemRPSD.h:86
void DrawAll() override
Definition: TotemRPSD.cc:261
void setTof(double e)
G4THitsCollection< TotemRPG4Hit > TotemRPG4HitCollection
Local3DPoint ConvertToLocal3DPoint(const G4ThreeVector &point) const
void setP(double e)
double Tof_
Definition: TotemRPSD.h:101
G4ThreeVector setToLocal(const G4ThreeVector &globalPoint)
Definition: TotemRPSD.cc:229
void createNewHit()
Definition: TotemRPSD.cc:189
void setIncidentEnergy(double e)
Definition: TotemRPG4Hit.cc:87
G4int tSliceID_
Definition: TotemRPSD.h:89
G4VPhysicalVolume * currentPV_
Definition: TotemRPSD.h:87
std::unique_ptr< TotemRPVDetectorOrganization > numberingScheme_
Definition: TotemRPSD.h:57
void storeHit(TotemRPG4Hit *)
Definition: TotemRPSD.cc:180
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
bool ProcessHits(G4Step *step, G4TouchableHistory *tHistory) override
Definition: TotemRPSD.cc:113
void PrintAll() override
Definition: TotemRPSD.cc:263
int unitID() const
Definition: TotemRPG4Hit.cc:92
int eventno_
Definition: TotemRPSD.h:113
void setVx(double p)
void setPy(double p)
TotemRPG4HitCollection * theHC_
Definition: TotemRPSD.h:83
double p() const
Definition: TotemRPG4Hit.cc:99
double Pabs_
Definition: TotemRPSD.h:99
G4double tSlice_
Definition: TotemRPSD.h:90
Geom::Theta< T > theta() const
Definition: PV3DBase.h:72
static constexpr double rp_garage_position_
Definition: TotemRPSD.h:44
void setLocalEntry(const G4ThreeVector &theLocalEntryPoint)
Definition: TotemRPG4Hit.cc:82
void stepInfo(const G4Step *aStep)
Definition: TotemRPSD.cc:124
void clearHits() override
Definition: TotemRPSD.cc:286
void printHitInfo()
Definition: TotemRPSD.cc:77
void Initialize(G4HCofThisEvent *HCE) override
Definition: TotemRPSD.cc:66
void EndOfEvent(G4HCofThisEvent *eventHC) override
Definition: TotemRPSD.cc:236
double phiAtEntry() const
std::string const collectionName[nCollections]
Definition: Collections.h:59
void setExit(G4ThreeVector xyz)
Definition: TotemRPG4Hit.cc:79
void clear() override
Definition: TotemRPSD.cc:259
G4ThreeVector theLocalEntryPoint_
Definition: TotemRPSD.h:96
int ParentId_
Definition: TotemRPSD.h:108
uint32_t setDetUnitId(const G4Step *step) override
Definition: TotemRPSD.cc:176
printf("params %d %f %f %f\n", minT, eps, errmax, chi2max)
short ParticleType_
Definition: TotemRPSD.h:103
G4StepPoint * postStepPoint_
Definition: TotemRPSD.h:93
int verbosity_
Definition: TotemRPSD.h:60
unsigned int unitID_
Definition: TotemRPSD.h:88
G4int primaryID_
Definition: TotemRPSD.h:89
void setLocalExit(const G4ThreeVector &theLocalExitPoint)
Definition: TotemRPG4Hit.cc:84
void setThetaAtEntry(double t)
void setNumberingScheme(TotemRPVDetectorOrganization *scheme)
Definition: TotemRPSD.cc:273
G4StepPoint * preStepPoint_
Definition: TotemRPSD.h:92
void setUnitID(unsigned int i)
Definition: TotemRPG4Hit.cc:93
double energyLoss() const
Log< level::Info, false > LogInfo
double PhiAtEntry_
Definition: TotemRPSD.h:106
double Vy_
Definition: TotemRPSD.h:109
G4ThreeVector theLocalExitPoint_
Definition: TotemRPSD.h:97
double ThetaAtEntry_
Definition: TotemRPSD.h:105
double Vz_
Definition: TotemRPSD.h:109
double thetaAtEntry() const
void setPhiAtEntry(double f)
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
bool isPrimary() const
double thePz_
Definition: TotemRPSD.h:100
static constexpr unsigned int maxTotemHits_
Definition: TotemRPSD.h:47
G4ThreeVector localExit() const
Definition: TotemRPG4Hit.cc:83
void setVz(double p)
double thePx_
Definition: TotemRPSD.h:100
double incidentEnergy_
Definition: TotemRPSD.h:79
double tof() const
void setParentId(int p)
TotemRPSD(const std::string &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
Definition: TotemRPSD.cc:33
void setTrackID(int i)
Definition: TotemRPG4Hit.cc:90
G4ThreeVector exitPoint_
Definition: TotemRPSD.h:95
void setPz(double p)
void update(const BeginOfEvent *) override
This routine will be called when the appropriate signal arrives.
Definition: TotemRPSD.cc:279
list entry
Definition: mps_splice.py:68
void summarize()
Definition: TotemRPSD.cc:257
void setParticleType(short i)
G4int hcID_
Definition: TotemRPSD.h:82
std::vector< PSimHit > PSimHitContainer
bool isPrimary(const G4Track *track)
Definition: TotemRPSD.cc:288
Log< level::Warning, false > LogWarning
~TotemRPSD() override
Definition: TotemRPSD.cc:64
double Eloss_
Definition: TotemRPSD.h:102
void setEnergyLoss(double e)
std::unique_ptr< TrackingSlaveSD > slave_
Definition: TotemRPSD.h:56
void setVy(double p)
double thePy_
Definition: TotemRPSD.h:100
void setPx(double p)
G4ThreeVector hitPoint_
Definition: TotemRPSD.h:94
TotemRPG4Hit * currentHit_
Definition: TotemRPSD.h:85
double Vx_
Definition: TotemRPSD.h:109
G4String name_
Definition: TotemRPSD.h:81
void setTimeSlice(double d)
Definition: TotemRPG4Hit.cc:96
#define LogDebug(id)