CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
MuonSensitiveDetector.cc
Go to the documentation of this file.
8 
10 
13 
17 
21 
22 #include "G4SDManager.hh"
23 #include "G4VProcess.hh"
24 #include "G4EventManager.hh"
25 #include "G4SystemOfUnits.hh"
26 
28 
29 #include <iostream>
30 
32  const DDCompactView & cpv,
33  const SensitiveDetectorCatalog & clg,
34  edm::ParameterSet const & p,
35  const SimTrackManager* manager)
36  : SensitiveTkDetector(name, cpv, clg, p),
37  thePV(0), theHit(0), theDetUnitId(0), theTrackID(0), theManager(manager)
38 {
39  edm::ParameterSet m_MuonSD = p.getParameter<edm::ParameterSet>("MuonSD");
40  STenergyPersistentCut = m_MuonSD.getParameter<double>("EnergyThresholdForPersistency");//Default 1. GeV
41  STallMuonsPersistent = m_MuonSD.getParameter<bool>("AllMuonsPersistent");
42  printHits = m_MuonSD.getParameter<bool>("PrintHits");
43 
44  //
45  // Here simply create 1 MuonSlaveSD for the moment
46  //
47 
48  LogDebug("MuonSimDebug") << "create MuonSubDetector "<<name<<std::endl;
49 
50  detector = new MuonSubDetector(name);
51 
52  LogDebug("MuonSimDebug") << "create MuonFrameRotation"<<std::endl;
53 
54  if (detector->isEndcap()) {
55  // cout << "MuonFrameRotation create MuonEndcapFrameRotation"<<endl;
57  } else if (detector->isRPC()) {
58  // cout << "MuonFrameRotation create MuonRPCFrameRotation"<<endl;
60  } else if (detector->isGEM()) {
61  // cout << "MuonFrameRotation create MuonGEMFrameRotation"<<endl;
63  } else if (detector->isME0()) {
64  // cout << "MuonFrameRotation create MuonME0FrameRotation"<<endl;
66  } else {
67  theRotation = 0;
68  }
69  LogDebug("MuonSimDebug") << "create MuonSlaveSD"<<std::endl;
71  LogDebug("MuonSimDebug") << "create MuonSimHitNumberingScheme"<<std::endl;
73  g4numbering = new MuonG4Numbering(cpv);
74 
75 
76  //
77  // Now attach the right detectors (LogicalVolumes) to me
78  //
79  const std::vector<std::string>& lvNames = clg.logicalNames(name);
80  this->Register();
81  for (std::vector<std::string>::const_iterator it = lvNames.begin(); it != lvNames.end(); it++){
82  LogDebug("MuonSimDebug") << name << " MuonSensitiveDetector:: attaching SD to LV " << *it << std::endl;
83  this->AssignSD(*it);
84  }
85 
86  if (printHits) {
87  thePrinter = new SimHitPrinter("HitPositionOSCAR.dat");
88  }
89 
90 
91  LogDebug("MuonSimDebug") << " EnergyThresholdForPersistency " << STenergyPersistentCut << " AllMuonsPersistent " << STallMuonsPersistent << std::endl;
92 
95 
96 }
97 
98 
100  delete g4numbering;
101  delete numbering;
102  delete slaveMuon;
103  delete theRotation;
104  delete detector;
105 
107 
108  delete myG4TrackToParticleID;
109 }
110 
112  clearHits();
113 
114  //----- Initialize variables to check if two steps belong to same hit
115  thePV = 0;
116  theDetUnitId = 0;
117  theTrackID = 0;
118 
119 }
120 
121 void MuonSensitiveDetector::update(const ::EndOfEvent * ev)
122 {
123  //slaveMuon->renumbering(theManager);
124 }
125 
126 
128 {
129  LogDebug("MuonSimDebug") << "MuonSensitiveDetector::clearHits"<<std::endl;
131 }
132 
133 bool MuonSensitiveDetector::ProcessHits(G4Step * aStep, G4TouchableHistory * ROhist)
134 {
135  LogDebug("MuonSimDebug") <<" MuonSensitiveDetector::ProcessHits "<<InitialStepPosition(aStep,WorldCoordinates)<<std::endl;
136 
137  // TimeMe t1( theHitTimer, false);
138 
139  if (aStep->GetTotalEnergyDeposit()>0.){
140  // do not count neutrals that are killed by User Limits MinEKine
141  if( aStep->GetTrack()->GetDynamicParticle()->GetCharge() != 0 ){
142 
143  if (newHit(aStep)) {
144  saveHit();
145  createHit(aStep);
146  } else {
147  updateHit(aStep);
148  }
149  return true;
150  } else {
151  storeVolumeAndTrack(aStep);
152  return false;
153  }
154  }
155  return false;
156 }
157 
158 uint32_t MuonSensitiveDetector::setDetUnitId(G4Step * aStep)
159 {
160  // G4VPhysicalVolume * pv = aStep->GetPreStepPoint()->GetPhysicalVolume();
162 
163  std::stringstream MuonBaseNumber;
164  // LogDebug :: Print MuonBaseNumber
165  MuonBaseNumber << "MuonNumbering :: number of levels = "<<num.getLevels()<<std::endl;
166  MuonBaseNumber << "Level \t SuperNo \t BaseNo"<<std::endl;
167  for (int level=1;level<=num.getLevels();level++) {
168  MuonBaseNumber << level << " \t " << num.getSuperNo(level)
169  << " \t " << num.getBaseNo(level) << std::endl;
170  }
171  std::string MuonBaseNumbr = MuonBaseNumber.str();
172 
173  LogDebug("MuonSimDebug") <<"MuonSensitiveDetector::setDetUnitId :: "<<MuonBaseNumbr;
174  LogDebug("MuonSimDebug") <<"MuonSensitiveDetector::setDetUnitId :: MuonDetUnitId = "<<(numbering->baseNumberToUnitNumber(num));
175  return numbering->baseNumberToUnitNumber(num);
176 }
177 
178 
180  if (theRotation !=0 ) {
181  return theRotation->transformPoint(in,s);
182  }
183  return (in);
184 }
185 
187  return Local3DPoint(in.x()/cm,in.y()/cm,in.z()/cm);
188 }
189 
191  return Global3DPoint(in.x()/cm,in.y()/cm,in.z()/cm);
192 }
193 
195  G4VPhysicalVolume* pv = aStep->GetPreStepPoint()->GetPhysicalVolume();
196  G4Track * t = aStep->GetTrack();
197  thePV=pv;
198  theTrackID=t->GetTrackID();
199 }
200 
201 bool MuonSensitiveDetector::newHit(G4Step * aStep){
202 
203  G4VPhysicalVolume* pv = aStep->GetPreStepPoint()->GetPhysicalVolume();
204  G4Track * t = aStep->GetTrack();
205  uint32_t currentUnitId=setDetUnitId(aStep);
206  LogDebug("MuonSimDebug") <<"MuonSensitiveDetector::newHit :: currentUnitId = "<<currentUnitId;
207  unsigned int currentTrackID=t->GetTrackID();
208  //unsigned int currentEventID=G4EventManager::GetEventManager()->GetConstCurrentEvent()->GetEventID();
209  bool changed=((pv!=thePV) ||
210  (currentUnitId!=theDetUnitId) ||
211  (currentTrackID!=theTrackID));
212  return changed;
213 }
214 
215 void MuonSensitiveDetector::createHit(G4Step * aStep){
216 
217  G4Track * theTrack = aStep->GetTrack();
218 
219  Local3DPoint theEntryPoint;
220  Local3DPoint theExitPoint;
221 
222  if (detector->isBarrel()) {
223  theEntryPoint= toOrcaUnits(toOrcaRef(InitialStepPositionVsParent(aStep,1),aStep)); // 1 level up
224  theExitPoint= toOrcaUnits(toOrcaRef(FinalStepPositionVsParent(aStep,1),aStep));
225  } else if (detector->isEndcap()) {
226  // save local z at current level
227  theEntryPoint= toOrcaUnits(toOrcaRef(InitialStepPosition(aStep,LocalCoordinates),aStep));
228  theExitPoint= toOrcaUnits(toOrcaRef(FinalStepPosition(aStep,LocalCoordinates),aStep));
229  float zentry = theEntryPoint.z();
230  float zexit = theExitPoint.z();
231  Local3DPoint tempEntryPoint= toOrcaUnits(toOrcaRef(InitialStepPositionVsParent(aStep,4),aStep)); // 4 levels up
232  Local3DPoint tempExitPoint= toOrcaUnits(toOrcaRef(FinalStepPositionVsParent(aStep,4),aStep));
233  // reset local z from z wrt deep-parent volume to z wrt low-level volume
234  theEntryPoint = Local3DPoint( tempEntryPoint.x(), tempEntryPoint.y(), zentry );
235  theExitPoint = Local3DPoint( tempExitPoint.x(), tempExitPoint.y(), zexit );
236  } else {
237  theEntryPoint= toOrcaUnits(toOrcaRef(InitialStepPosition(aStep,LocalCoordinates),aStep));
238  theExitPoint= toOrcaUnits(toOrcaRef(FinalStepPosition(aStep,LocalCoordinates),aStep));
239  }
240 
241  float thePabs = aStep->GetPreStepPoint()->GetMomentum().mag()/GeV;
242  float theTof = aStep->GetPreStepPoint()->GetGlobalTime()/nanosecond;
243  float theEnergyLoss = aStep->GetTotalEnergyDeposit()/GeV;
244  // int theParticleType = theTrack->GetDefinition()->GetPDGEncoding();
245  int theParticleType = myG4TrackToParticleID->particleID(theTrack);
246  G4ThreeVector gmd = aStep->GetPreStepPoint()->GetMomentumDirection();
247  G4ThreeVector lmd = ((const G4TouchableHistory *)(aStep->GetPreStepPoint()->GetTouchable()))->GetHistory()
248  ->GetTopTransform().TransformAxis(gmd);
249  Local3DPoint lnmd = toOrcaRef(ConvertToLocal3DPoint(lmd),aStep);
250  float theThetaAtEntry = lnmd.theta();
251  float thePhiAtEntry = lnmd.phi();
252 
253  storeVolumeAndTrack( aStep );
254  theDetUnitId = setDetUnitId(aStep);
255 
256  Global3DPoint theGlobalPos;
257  if (printHits) {
258  Local3DPoint theGlobalHelp = InitialStepPosition(aStep,WorldCoordinates);
259  theGlobalEntry = toOrcaUnits(Global3DPoint (theGlobalHelp.x(),theGlobalHelp.y(),theGlobalHelp.z()));
260 
261  G4StepPoint * preStepPoint = aStep->GetPreStepPoint();
262  const G4TouchableHistory * theTouchable=(const G4TouchableHistory *)
263  (preStepPoint->GetTouchable());
264  theGlobalHelp=ConvertToLocal3DPoint(theTouchable->GetTranslation());
265  theGlobalPos = toOrcaUnits(Global3DPoint (theGlobalHelp.x(),theGlobalHelp.y(),theGlobalHelp.z()));
266  // const G4RotationMatrix * theGlobalRot = theTouchable->GetRotation();
267  }
268 
269  LogDebug("MuonSimDebug") << "MuonSensitiveDetector::createHit UpdatablePSimHit"<<std::endl;
270 
271  theHit = new UpdatablePSimHit(theEntryPoint,theExitPoint,thePabs,theTof,
272  theEnergyLoss,theParticleType,theDetUnitId,
273  theTrackID,theThetaAtEntry,thePhiAtEntry,
274  theG4ProcessTypeEnumerator->processId(theTrack->GetCreatorProcess()));
275 
276 
277  LogDebug("MuonSimDebug") <<"=== NEW ==================> ELOSS = "<<theEnergyLoss<<" "
278  <<thePV->GetLogicalVolume()->GetName()<<std::endl;
279  const G4VProcess* p = aStep->GetPostStepPoint()->GetProcessDefinedStep();
280  const G4VProcess* p2 = aStep->GetPreStepPoint()->GetProcessDefinedStep();
281  if (p)
282  LogDebug("MuonSimDebug") <<" POST PROCESS = "<<p->GetProcessName()<<std::endl;
283  if (p2)
284  LogDebug("MuonSimDebug") <<" PRE PROCESS = "<<p2->GetProcessName()<<std::endl;
285  LogDebug("MuonSimDebug") << "newhit theta " << theThetaAtEntry<<std::endl;
286  LogDebug("MuonSimDebug") << "newhit phi " << thePhiAtEntry<<std::endl;
287  LogDebug("MuonSimDebug") << "newhit pabs " << thePabs<<std::endl;
288  LogDebug("MuonSimDebug") << "newhit tof " << theTof<<std::endl;
289  LogDebug("MuonSimDebug") << "newhit track " << theTrackID<<std::endl;
290  LogDebug("MuonSimDebug") << "newhit entry " << theEntryPoint<<std::endl;
291  LogDebug("MuonSimDebug") << "newhit exit " << theExitPoint<<std::endl;
292  LogDebug("MuonSimDebug") << "newhit eloss " << theEnergyLoss << std::endl;
293  LogDebug("MuonSimDebug") << "newhit detid " << theDetUnitId<<std::endl;
294  LogDebug("MuonSimDebug") << "newhit delta " << (theExitPoint-theEntryPoint)<<std::endl;
295  LogDebug("MuonSimDebug") << "newhit deltu " << (theExitPoint-theEntryPoint).unit();
296  LogDebug("MuonSimDebug") << " " << (theExitPoint-theEntryPoint).mag()<<std::endl;
297  LogDebug("MuonSimDebug") << "newhit glob " << theGlobalEntry<<std::endl;
298  LogDebug("MuonSimDebug") << "newhit dpos " << theGlobalPos<<std::endl;
299  LogDebug("MuonSimDebug") << "newhit drot " << std::endl;
300  // theGlobalRot->print(LogDebug("MuonSimDebug"));
301 
302 
303  //
304  //----- SimTracks: Make it persistent?
305  //
306  int thePID = theTrack->GetDefinition()->GetPDGEncoding();
307  LogDebug("MuonSimDebug") << " checking simtrack " << thePID << " " << thePabs << " STenergyPersistentCut " << STenergyPersistentCut << std::endl;
308 
309  if( thePabs*GeV > STenergyPersistentCut
310  || ( abs(thePID) == 13 && STallMuonsPersistent ) ){
312  LogDebug("MuonSimDebug") <<" track leaving hit in muons made selected for persistency"<<std::endl;
313 
314  info->storeTrack(true);
315  }
316 
317 }
318 
319 void MuonSensitiveDetector::updateHit(G4Step * aStep){
320  // float thePabs = aStep->GetPreStepPoint()->GetMomentum().mag()/GeV;
321  // Local3DPoint theEntryPoint= InitialStepPosition(aStep,LocalCoordinates);
322 
323 
324  Local3DPoint theExitPoint;
325 
326  if (detector->isBarrel()) {
327  theExitPoint= toOrcaUnits(toOrcaRef(FinalStepPositionVsParent(aStep,1),aStep));
328  } else if (detector->isEndcap()) {
329  // save local z at current level
330  theExitPoint= toOrcaUnits(toOrcaRef(FinalStepPosition(aStep,LocalCoordinates),aStep));
331  float zexit = theExitPoint.z();
332  Local3DPoint tempExitPoint= toOrcaUnits(toOrcaRef(FinalStepPositionVsParent(aStep,4),aStep));
333  theExitPoint = Local3DPoint( tempExitPoint.x(), tempExitPoint.y(), zexit );
334  } else {
335  theExitPoint= toOrcaUnits(toOrcaRef(FinalStepPosition(aStep,LocalCoordinates),aStep));
336  }
337 
338  float theEnergyLoss = aStep->GetTotalEnergyDeposit()/GeV;
339 
340  if( theHit == 0 ){
341  std::cerr << "!!ERRROR in MuonSensitiveDetector::updateHit. It is called when there is no hit " << std::endl;
342  }
343 
344  theHit->updateExitPoint(theExitPoint);
345  theHit->addEnergyLoss(theEnergyLoss);
346 
347  LogDebug("MuonSimDebug") <<"=== UPDATE ===============> ELOSS = "<<theEnergyLoss<<" "
348  <<thePV->GetLogicalVolume()->GetName()<<std::endl;
349  const G4VProcess* p = aStep->GetPostStepPoint()->GetProcessDefinedStep();
350  const G4VProcess* p2 = aStep->GetPreStepPoint()->GetProcessDefinedStep();
351  if (p)
352  LogDebug("MuonSimDebug") <<" POST PROCESS = "<<p->GetProcessName()<<std::endl;
353  if (p2)
354  LogDebug("MuonSimDebug") <<" PRE PROCESS = "<<p2->GetProcessName()<<std::endl;
355  LogDebug("MuonSimDebug") << "updhit exit " << theExitPoint<<std::endl;
356  LogDebug("MuonSimDebug") << "updhit eloss " << theHit->energyLoss() <<std::endl;
357  LogDebug("MuonSimDebug") << "updhit detid " << theDetUnitId<<std::endl;
358  LogDebug("MuonSimDebug") << "updhit delta " << (theExitPoint-theHit->entryPoint())<<std::endl;
359  LogDebug("MuonSimDebug") << "updhit deltu " << (theExitPoint-theHit->entryPoint()).unit();
360  LogDebug("MuonSimDebug") << " " << (theExitPoint-theHit->entryPoint()).mag()<<std::endl;
361 
362 }
363 
365 
366  if (theHit) {
367  if (printHits) {
370  // thePrinter->printTrack(theHit->trackId());
371  //thePrinter->printPabs(theHit->pabs());
372  //thePrinter->printEloss(theHit->energyLoss());
375  }
377  // seems the hit does not want to be deleted
378  // done by the hit collection?
379  delete theHit;
380  theHit = 0; //set it to 0, because you are checking that is 0
381  }
382 
383 }
384 
386 {
387  G4VUserTrackInformation* temp = gTrack->GetUserInformation();
388  if (temp == 0){
389  std::cerr <<" ERROR: no G4VUserTrackInformation available"<<std::endl;
390  abort();
391  }else{
392  TrackInformation* info = dynamic_cast<TrackInformation*>(temp);
393  if (info ==0){
394  std::cerr <<" ERROR: TkSimTrackSelection: the UserInformation does not appear to be a TrackInformation"<<std::endl;
395  abort();
396  }
397  return info;
398  }
399 }
400 
401 void MuonSensitiveDetector::EndOfEvent(G4HCofThisEvent*)
402 {
403 // TimeMe t("MuonSensitiveDetector::EndOfEvent", false);
404  // LogDebug("MuonSimDebug") << "MuonSensitiveDetector::EndOfEvent saving last hit en event " << std::endl;
405  saveHit();
406 }
407 
408 
410  //
411  // do it once for low, once for High
412  //
413 
414  if (slaveMuon->name() == n) c=slaveMuon->hits();
415 
416 }
417 
418 std::vector<std::string> MuonSensitiveDetector::getNames(){
419  std::vector<std::string> temp;
420  temp.push_back(slaveMuon->name());
421  return temp;
422 }
423 
425 
426  G4StepPoint * preStepPoint = currentStep->GetPreStepPoint();
427  G4ThreeVector globalCoordinates = preStepPoint->GetPosition();
428 
429  const G4TouchableHistory * theTouchable=(const G4TouchableHistory *)
430  (preStepPoint->GetTouchable());
431 
432  G4int depth = theTouchable->GetHistory()->GetDepth();
433  G4ThreeVector localCoordinates = theTouchable->GetHistory()
434  ->GetTransform(depth-levelsUp).TransformPoint(globalCoordinates);
435 
436  return ConvertToLocal3DPoint(localCoordinates);
437 }
438 
439 Local3DPoint MuonSensitiveDetector::FinalStepPositionVsParent(G4Step * currentStep, G4int levelsUp) {
440 
441  G4StepPoint * postStepPoint = currentStep->GetPostStepPoint();
442  G4StepPoint * preStepPoint = currentStep->GetPreStepPoint();
443  G4ThreeVector globalCoordinates = postStepPoint->GetPosition();
444 
445  const G4TouchableHistory * theTouchable = (const G4TouchableHistory *)
446  (preStepPoint->GetTouchable());
447 
448  G4int depth = theTouchable->GetHistory()->GetDepth();
449  G4ThreeVector localCoordinates = theTouchable->GetHistory()
450  ->GetTransform(depth-levelsUp).TransformPoint(globalCoordinates);
451 
452  return ConvertToLocal3DPoint(localCoordinates);
453 }
#define LogDebug(id)
T getParameter(std::string const &) const
tuple t
Definition: tree.py:139
int i
Definition: DBlmapReader.cc:9
Local3DPoint FinalStepPositionVsParent(G4Step *currentStep, G4int levelsUp)
static const TGPicture * info(bool iBackgroundIsBlack)
const double GeV
Definition: MathUtil.h:16
bool storeTrack() const
void printLocal(LocalPoint, LocalPoint) const
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
std::string name() const
int getBaseNo(int level) const
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
T y() const
Definition: PV3DBase.h:63
bool ev
void fillHits(edm::PSimHitContainer &, std::string use)
virtual Local3DPoint transformPoint(const Local3DPoint &, const G4Step *) const =0
std::string name()
const std::vector< std::string > & logicalNames(const std::string &readoutName) const
type of data representation of DDCompactView
Definition: DDCompactView.h:76
virtual void AssignSD(const std::string &vname)
MuonSubDetector * detector
virtual void EndOfEvent(G4HCofThisEvent *)
unsigned int processId(const G4VProcess *p)
Geom::Theta< T > theta() const
Definition: PV3DBase.h:75
std::vector< PSimHit > & hits()
Point3DBase< float, GlobalTag > Global3DPoint
Definition: GlobalPoint.h:7
Local3DPoint exitPoint() const
Exit point in the local Det frame.
Definition: PSimHit.h:38
Local3DPoint ConvertToLocal3DPoint(const G4ThreeVector &point)
void updateExitPoint(const Local3DPoint &exit)
virtual void Initialize()
G4ProcessTypeEnumerator * theG4ProcessTypeEnumerator
string unit
Definition: csvLumiCalc.py:46
T z() const
Definition: PV3DBase.h:64
Local3DPoint toOrcaUnits(Local3DPoint)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< std::string > getNames()
double p2[4]
Definition: TauolaWrapper.h:90
const SimTrackManager * theManager
void printGlobal(GlobalPoint) const
static int particleID(const G4Track *)
MuonFrameRotation * theRotation
virtual int baseNumberToUnitNumber(const MuonBaseNumber &)
MuonBaseNumber PhysicalVolumeToBaseNumber(const G4Step *aStep)
Point3DBase< float, LocalTag > Local3DPoint
Definition: LocalPoint.h:9
void addEnergyLoss(float eloss)
void startNewSimHit(std::string)
MuonSimHitNumberingScheme * numbering
Local3DPoint toOrcaRef(Local3DPoint in, G4Step *s)
G4TrackToParticleID * myG4TrackToParticleID
void printId(int) const
int getLevels() const
float energyLoss() const
The energy deposit in the PSimHit, in ???.
Definition: PSimHit.h:75
MuonSensitiveDetector(std::string, const DDCompactView &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
int getSuperNo(int level) const
G4VPhysicalVolume * thePV
virtual bool processHits(const PSimHit &)
void update(const BeginOfEvent *)
This routine will be called when the appropriate signal arrives.
std::vector< PSimHit > PSimHitContainer
tuple level
Definition: testEve_cfg.py:34
Local3DPoint InitialStepPositionVsParent(G4Step *currentStep, G4int levelsUp)
Local3DPoint FinalStepPosition(G4Step *s, coordinates)
virtual uint32_t setDetUnitId(G4Step *)
T x() const
Definition: PV3DBase.h:62
TrackInformation * getOrCreateTrackInformation(const G4Track *theTrack)
Local3DPoint entryPoint() const
Entry point in the local Det frame.
Definition: PSimHit.h:35
Local3DPoint InitialStepPosition(G4Step *s, coordinates)
unsigned int detUnitId() const
Definition: PSimHit.h:93
MuonG4Numbering * g4numbering
virtual G4bool ProcessHits(G4Step *, G4TouchableHistory *)