CMS 3D CMS Logo

MuonSensitiveDetector.cc
Go to the documentation of this file.
9 
13 
17 
21 
22 #include "G4VProcess.hh"
23 #include "G4EventManager.hh"
24 #include "G4SystemOfUnits.hh"
25 #include "G4Step.hh"
26 #include "G4StepPoint.hh"
27 #include "G4Track.hh"
28 
30 
31 #include <iostream>
32 
33 //#define DebugLog
34 
36  const DDCompactView& cpv,
37  const SensitiveDetectorCatalog& clg,
38  edm::ParameterSet const& p,
39  const SimTrackManager* manager)
40  : SensitiveTkDetector(name, cpv, clg, p),
41  thePV(nullptr),
42  theHit(nullptr),
43  theDetUnitId(0),
44  newDetUnitId(0),
45  theTrackID(0),
46  theManager(manager) {
47  edm::ParameterSet m_MuonSD = p.getParameter<edm::ParameterSet>("MuonSD");
48  ePersistentCutGeV = m_MuonSD.getParameter<double>("EnergyThresholdForPersistency") / CLHEP::GeV; //Default 1. GeV
49  allMuonsPersistent = m_MuonSD.getParameter<bool>("AllMuonsPersistent");
50  printHits = m_MuonSD.getParameter<bool>("PrintHits");
51 
52  //
53  // Here simply create 1 MuonSlaveSD for the moment
54  //
55  LogDebug("MuonSimDebug") << "create MuonSubDetector " << name;
56  detector = new MuonSubDetector(name);
57 
58  //The constants take time to calculate and are needed by many helpers
60  G4String sdet = "unknown";
61  if (detector->isEndcap()) {
63  sdet = "Endcap";
64  } else if (detector->isRPC()) {
65  theRotation = new MuonRPCFrameRotation(constants);
66  sdet = "RPC";
67  } else if (detector->isGEM()) {
68  theRotation = new MuonGEMFrameRotation(constants);
69  sdet = "GEM";
70  } else if (detector->isME0()) {
71  theRotation = new MuonME0FrameRotation(constants);
72  sdet = "ME0";
73  } else {
75  }
78  g4numbering = new MuonG4Numbering(constants);
79 
80  if (printHits) {
81  thePrinter = new SimHitPrinter("HitPositionOSCAR.dat");
82  }
83 
84  edm::LogVerbatim("MuonSensitiveDetector")
85  << " of type " << sdet << " <" << GetName() << "> EnergyThresholdForPersistency(GeV) "
86  << ePersistentCutGeV / CLHEP::GeV << " allMuonsPersistent: " << allMuonsPersistent;
87 
89 }
90 
92  delete g4numbering;
93  delete numbering;
94  delete slaveMuon;
95  delete theRotation;
96  delete detector;
98 }
99 
101  clearHits();
102  //----- Initialize variables to check if two steps belong to same hit
103  thePV = nullptr;
104  theDetUnitId = 0;
105  theTrackID = 0;
106 }
107 
109  LogDebug("MuonSimDebug") << "MuonSensitiveDetector::clearHits";
111 }
112 
113 bool MuonSensitiveDetector::ProcessHits(G4Step* aStep, G4TouchableHistory* ROhist) {
114  LogDebug("MuonSimDebug") << " MuonSensitiveDetector::ProcessHits " << InitialStepPosition(aStep, WorldCoordinates);
115 
116  if (aStep->GetTotalEnergyDeposit() > 0.) {
117  newDetUnitId = setDetUnitId(aStep);
118 
119  // do not count neutrals that are killed by User Limits MinEKine
120  //---VI: this is incorrect, neutral particle, like neutron may have local
121  // energy deposit, which potentially may make a hit
122  if (aStep->GetTrack()->GetDynamicParticle()->GetCharge() != 0) {
123  if (newHit(aStep)) {
124  saveHit();
125  createHit(aStep);
126  } else {
127  updateHit(aStep);
128  }
129  } else {
130  thePV = aStep->GetPreStepPoint()->GetPhysicalVolume();
131  theTrackID = aStep->GetTrack()->GetTrackID();
133  }
134  }
135  return true;
136 }
137 
138 uint32_t MuonSensitiveDetector::setDetUnitId(const G4Step* aStep) {
140 
141 #ifdef DebugLog
142  std::stringstream MuonBaseNumber;
143  MuonBaseNumber << "MuonNumbering :: number of levels = " << num.getLevels() << std::endl;
144  MuonBaseNumber << "Level \t SuperNo \t BaseNo" << std::endl;
145  for (int level = 1; level <= num.getLevels(); level++) {
146  MuonBaseNumber << level << " \t " << num.getSuperNo(level) << " \t " << num.getBaseNo(level) << std::endl;
147  }
148  std::string MuonBaseNumbr = MuonBaseNumber.str();
149 
150  LogDebug("MuonSimDebug") << "MuonSensitiveDetector::setDetUnitId :: " << MuonBaseNumbr;
151  LogDebug("MuonSimDebug") << "MuonSensitiveDetector::setDetUnitId :: MuonDetUnitId = "
153 #endif
154  return numbering->baseNumberToUnitNumber(num);
155 }
156 
157 bool MuonSensitiveDetector::newHit(const G4Step* aStep) {
158  return (!theHit || (aStep->GetTrack()->GetTrackID() != theTrackID) ||
159  (aStep->GetPreStepPoint()->GetPhysicalVolume() != thePV) || newDetUnitId != theDetUnitId);
160 }
161 
162 void MuonSensitiveDetector::createHit(const G4Step* aStep) {
163  Local3DPoint theEntryPoint;
164  Local3DPoint theExitPoint;
165 
166  if (detector->isBarrel()) {
167  // 1 levels up
168  theEntryPoint = cmsUnits(theRotation->transformPoint(InitialStepPositionVsParent(aStep, 1), aStep));
169  theExitPoint = cmsUnits(theRotation->transformPoint(FinalStepPositionVsParent(aStep, 1), aStep));
170  } else if (detector->isEndcap()) {
171  // save local z at current level
172  theEntryPoint = theRotation->transformPoint(InitialStepPosition(aStep, LocalCoordinates), aStep);
173  theExitPoint = theRotation->transformPoint(FinalStepPosition(aStep, LocalCoordinates), aStep);
174  float zentry = theEntryPoint.z();
175  float zexit = theExitPoint.z();
176  // 4 levels up
179  // reset local z from z wrt deep-parent volume to z wrt low-level volume
180  theEntryPoint = cmsUnits(Local3DPoint(tempEntry.x(), tempEntry.y(), zentry));
181  theExitPoint = cmsUnits(Local3DPoint(tempExit.x(), tempExit.y(), zexit));
182  } else {
183  theEntryPoint = cmsUnits(theRotation->transformPoint(InitialStepPosition(aStep, LocalCoordinates), aStep));
184  theExitPoint = cmsUnits(theRotation->transformPoint(FinalStepPosition(aStep, LocalCoordinates), aStep));
185  }
186 
187  const G4Track* theTrack = aStep->GetTrack();
188  const G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
189 
190  float thePabs = preStepPoint->GetMomentum().mag() / CLHEP::GeV;
191  float theTof = preStepPoint->GetGlobalTime() / CLHEP::nanosecond;
192  float theEnergyLoss = aStep->GetTotalEnergyDeposit() / CLHEP::GeV;
193  int theParticleType = G4TrackToParticleID::particleID(theTrack);
194 
196  thePV = preStepPoint->GetPhysicalVolume();
197  theTrackID = theTrack->GetTrackID();
198 
199  // convert momentum direction it to local frame
200  const G4ThreeVector& gmd = preStepPoint->GetMomentumDirection();
201  G4ThreeVector lmd = static_cast<const G4TouchableHistory*>(preStepPoint->GetTouchable())
202  ->GetHistory()
203  ->GetTopTransform()
204  .TransformAxis(gmd);
206  lnmd = theRotation->transformPoint(lnmd, aStep);
207  float theThetaAtEntry = lnmd.theta();
208  float thePhiAtEntry = lnmd.phi();
209 
210  theHit = new UpdatablePSimHit(theEntryPoint,
211  theExitPoint,
212  thePabs,
213  theTof,
214  theEnergyLoss,
215  theParticleType,
216  theDetUnitId,
217  theTrackID,
218  theThetaAtEntry,
219  thePhiAtEntry,
220  theG4ProcessTypeEnumerator->processId(theTrack->GetCreatorProcess()));
221 
222  // Make track persistent
223  int thePID = std::abs(theTrack->GetDefinition()->GetPDGEncoding());
224  //---VI - in parameters cut in energy is declared but applied to momentum
225  if (thePabs > ePersistentCutGeV || (thePID == 13 && allMuonsPersistent)) {
227  info->storeTrack(true);
228  }
229 
230 #ifdef DebugLog
231  edm::LogVerbatim("MuonSimDebug") << "=== NEW Muon hit for " << GetName() << " Edep(GeV)= " << theEnergyLoss << " "
232  << thePV->GetLogicalVolume()->GetName();
233  const G4VProcess* p = aStep->GetPostStepPoint()->GetProcessDefinedStep();
234  const G4VProcess* p2 = preStepPoint->GetProcessDefinedStep();
235  G4String sss = "";
236  if (p)
237  sss += " POST PROCESS: " + p->GetProcessName();
238  if (p2)
239  sss += "; PRE PROCESS: " + p2->GetProcessName();
240  if ("" != sss)
241  edm::LogVerbatim("MuonSimDebug") << sss;
242  edm::LogVerbatim("MuonSimDebug") << " theta= " << theThetaAtEntry << " phi= " << thePhiAtEntry << " Pabs(GeV/c) "
243  << thePabs << " Eloss(GeV)= " << theEnergyLoss << " Tof(ns)= " << theTof
244  << " trackID= " << theTrackID << " detID= " << theDetUnitId << "\n Local: entry "
245  << theEntryPoint << " exit " << theExitPoint << " delta "
246  << (theExitPoint - theEntryPoint) << "\n Global: entry "
247  << preStepPoint->GetPosition() << " exit "
248  << aStep->GetPostStepPoint()->GetPosition();
249 #endif
250 }
251 
252 void MuonSensitiveDetector::updateHit(const G4Step* aStep) {
253  Local3DPoint theExitPoint;
254 
255  if (detector->isBarrel()) {
256  theExitPoint = cmsUnits(theRotation->transformPoint(FinalStepPositionVsParent(aStep, 1), aStep));
257  } else if (detector->isEndcap()) {
258  // save local z at current level
259  theExitPoint = theRotation->transformPoint(FinalStepPosition(aStep, LocalCoordinates), aStep);
260  float zexit = theExitPoint.z();
261  Local3DPoint tempExitPoint = theRotation->transformPoint(FinalStepPositionVsParent(aStep, 4), aStep);
262  theExitPoint = cmsUnits(Local3DPoint(tempExitPoint.x(), tempExitPoint.y(), zexit));
263  } else {
264  theExitPoint = cmsUnits(theRotation->transformPoint(FinalStepPosition(aStep, LocalCoordinates), aStep));
265  }
266 
267  float theEnergyLoss = aStep->GetTotalEnergyDeposit() / CLHEP::GeV;
268 
269  theHit->updateExitPoint(theExitPoint);
270  theHit->addEnergyLoss(theEnergyLoss);
271 
272 #ifdef DebugLog
273  edm::LogVerbatim("MuonSimDebug") << "=== NEW Update muon hit for " << GetName() << " Edep(GeV)= " << theEnergyLoss
274  << " " << thePV->GetLogicalVolume()->GetName();
275  const G4VProcess* p = aStep->GetPostStepPoint()->GetProcessDefinedStep();
276  const G4VProcess* p2 = preStepPoint->GetProcessDefinedStep();
277  G4String sss = "";
278  if (p)
279  sss += " POST PROCESS: " + p->GetProcessName();
280  if (p2)
281  sss += "; PRE PROCESS: " + p2->GetProcessName();
282  if ("" != sss)
283  edm::LogVerbatim("MuonSimDebug") << sss;
284  edm::LogVerbatim("MuonSimDebug") << " delEloss(GeV)= " << theEnergyLoss << " Tof(ns)= " << theTof
285  << " trackID= " << theTrackID << " detID= " << theDetUnitId << " exit "
286  << theExitPoint;
287 #endif
288 }
289 
291  if (theHit) {
292  if (printHits) {
296  }
297  // hit is included into hit collection
299  delete theHit;
300  theHit = nullptr;
301  }
302 }
303 
304 void MuonSensitiveDetector::EndOfEvent(G4HCofThisEvent*) { saveHit(); }
305 
307  if (slaveMuon->name() == hname) {
308  cc = slaveMuon->hits();
309  }
310 }
311 
312 Local3DPoint MuonSensitiveDetector::InitialStepPositionVsParent(const G4Step* currentStep, G4int levelsUp) {
313  const G4StepPoint* preStepPoint = currentStep->GetPreStepPoint();
314  const G4ThreeVector& globalCoordinates = preStepPoint->GetPosition();
315 
316  const G4TouchableHistory* theTouchable = (const G4TouchableHistory*)(preStepPoint->GetTouchable());
317 
318  G4int depth = theTouchable->GetHistory()->GetDepth();
319  G4ThreeVector localCoordinates =
320  theTouchable->GetHistory()->GetTransform(depth - levelsUp).TransformPoint(globalCoordinates);
321 
322  return ConvertToLocal3DPoint(localCoordinates);
323 }
324 
325 Local3DPoint MuonSensitiveDetector::FinalStepPositionVsParent(const G4Step* currentStep, G4int levelsUp) {
326  const G4StepPoint* postStepPoint = currentStep->GetPostStepPoint();
327  const G4StepPoint* preStepPoint = currentStep->GetPreStepPoint();
328  const G4ThreeVector& globalCoordinates = postStepPoint->GetPosition();
329 
330  const G4TouchableHistory* theTouchable = (const G4TouchableHistory*)(preStepPoint->GetTouchable());
331 
332  G4int depth = theTouchable->GetHistory()->GetDepth();
333  G4ThreeVector localCoordinates =
334  theTouchable->GetHistory()->GetTransform(depth - levelsUp).TransformPoint(globalCoordinates);
335 
336  return ConvertToLocal3DPoint(localCoordinates);
337 }
#define LogDebug(id)
T getParameter(std::string const &) const
const G4VPhysicalVolume * thePV
static const TGPicture * info(bool iBackgroundIsBlack)
const double GeV
Definition: MathUtil.h:16
Local3DPoint FinalStepPositionVsParent(const G4Step *currentStep, G4int levelsUp)
bool storeTrack() const
G4bool ProcessHits(G4Step *, G4TouchableHistory *) override
Local3DPoint ConvertToLocal3DPoint(const G4ThreeVector &point) const
void updateHit(const G4Step *)
void printLocal(LocalPoint, LocalPoint) const
std::string name() const
int getBaseNo(int level) const
#define nullptr
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
T y() const
Definition: PV3DBase.h:63
std::string name()
Compact representation of the geometrical detector hierarchy.
Definition: DDCompactView.h:80
MuonSubDetector * detector
void fillHits(edm::PSimHitContainer &, const std::string &) override
Geom::Theta< T > theta() const
Definition: PV3DBase.h:75
Local3DPoint cmsUnits(const Local3DPoint &v)
std::vector< PSimHit > & hits()
unsigned int processId(const G4VProcess *p) const
Local3DPoint exitPoint() const
Exit point in the local Det frame.
Definition: PSimHit.h:46
T mag() const
Definition: PV3DBase.h:67
void createHit(const G4Step *)
void updateExitPoint(const Local3DPoint &exit)
virtual void Initialize()
G4ProcessTypeEnumerator * theG4ProcessTypeEnumerator
Local3DPoint FinalStepPosition(const G4Step *step, coordinates) const
T z() const
Definition: PV3DBase.h:64
bool newHit(const G4Step *)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
virtual Local3DPoint transformPoint(const Local3DPoint &, const G4Step *) const
double p2[4]
Definition: TauolaWrapper.h:90
const SimTrackManager * theManager
TrackInformation * cmsTrackInformation(const G4Track *aTrack)
static int particleID(const G4Track *)
MuonFrameRotation * theRotation
MuonSensitiveDetector(const std::string &, const DDCompactView &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
MuonBaseNumber PhysicalVolumeToBaseNumber(const G4Step *aStep)
uint32_t setDetUnitId(const G4Step *) override
Point3DBase< float, LocalTag > Local3DPoint
Definition: LocalPoint.h:9
void addEnergyLoss(float eloss)
void startNewSimHit(std::string)
void EndOfEvent(G4HCofThisEvent *) override
MuonSimHitNumberingScheme * numbering
void printId(int) const
int getLevels() const
void update(const BeginOfEvent *) override
This routine will be called when the appropriate signal arrives.
int getSuperNo(int level) const
virtual bool processHits(const PSimHit &)
std::vector< PSimHit > PSimHitContainer
T x() const
Definition: PV3DBase.h:62
int baseNumberToUnitNumber(const MuonBaseNumber &) override
Local3DPoint entryPoint() const
Entry point in the local Det frame.
Definition: PSimHit.h:43
unsigned int detUnitId() const
Definition: PSimHit.h:97
Local3DPoint InitialStepPositionVsParent(const G4Step *currentStep, G4int levelsUp)
UpdatablePSimHit * theHit
MuonG4Numbering * g4numbering
Local3DPoint InitialStepPosition(const G4Step *step, coordinates) const