CMS 3D CMS Logo

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