CMS 3D CMS Logo

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