CMS 3D CMS Logo

Bcm1fSD.cc
Go to the documentation of this file.
2 
5 
9 
12 
14 
19 
20 #include "G4Track.hh"
21 #include "G4SDManager.hh"
22 #include "G4VProcess.hh"
23 #include "G4EventManager.hh"
24 #include "G4Event.hh"
25 #include "G4VProcess.hh"
26 
27 #include <string>
28 #include <vector>
29 #include <iostream>
30 
31 #include "CLHEP/Units/GlobalSystemOfUnits.h"
32 
34  const SensitiveDetectorCatalog & clg,
35  edm::ParameterSet const & p,
36  const SimTrackManager* manager) :
37  SensitiveTkDetector(name, cpv, clg, p), mySimHit(nullptr),
38  oldVolume(nullptr), lastId(0), lastTrack(0), eventno(0) {
39 
40  edm::ParameterSet m_TrackerSD = p.getParameter<edm::ParameterSet>("Bcm1fSD");
41  energyCut = m_TrackerSD.getParameter<double>("EnergyThresholdForPersistencyInGeV")*GeV; //default must be 0.5 (?)
42  energyHistoryCut = m_TrackerSD.getParameter<double>("EnergyThresholdForHistoryInGeV")*GeV;//default must be 0.05 (?)
43 
44  edm::LogInfo("Bcm1fSD") <<" Criteria for Saving Tracker SimTracks: \n "
45  <<" History: "<<energyHistoryCut<< " MeV ; Persistency: "
46  << energyCut<<" MeV\n" <<" Constructing a Bcm1fSD";
47 
48  slave = new TrackingSlaveSD(name);
49 
52 }
53 
55  delete slave;
57  delete myG4TrackToParticleID;
58 }
59 
60 bool Bcm1fSD::ProcessHits(G4Step * aStep, G4TouchableHistory *) {
61 
62  LogDebug("Bcm1fSD") << " Entering a new Step "
63  << aStep->GetTotalEnergyDeposit() << " "
64  << aStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetName();
65 
66  G4Track * gTrack = aStep->GetTrack();
67  if ((unsigned int)(gTrack->GetTrackID()) != lastTrack) {
68 
69  if (gTrack->GetKineticEnergy() > energyCut){
71  info->storeTrack(true);
72  }
73  if (gTrack->GetKineticEnergy() > energyHistoryCut){
75  info->putInHistory();
76  }
77  }
78 
79  if (aStep->GetTotalEnergyDeposit()>0.) {
80  if (newHit(aStep) == true) {
81  sendHit();
82  createHit(aStep);
83  } else {
84  updateHit(aStep);
85  }
86  return true;
87  }
88  return false;
89 }
90 
91 uint32_t Bcm1fSD::setDetUnitId(const G4Step * aStep ) {
92 
93  unsigned int detId = 0;
94 
95  LogDebug("Bcm1fSD")<< " DetID = "<<detId;
96 
97  //Find number of levels
98  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
99  int level = 0;
100  if (touch) level = ((touch->GetHistoryDepth())+1);
101 
102  //Get name and copy numbers
103  if ( level > 1 ) {
104 
105  G4String sensorName = touch->GetVolume(0)->GetName();
106  G4String diamondName = touch->GetVolume(1)->GetName();
107  G4String detectorName = touch->GetVolume(2)->GetName();
108  G4String volumeName = touch->GetVolume(3)->GetName();
109 
110  if ( sensorName != "BCM1FSensor" )
111  std::cout << " Bcm1fSD::setDetUnitId -w- Sensor name not BCM1FSensor " << std::endl;
112  if ( detectorName != "BCM1F" )
113  std::cout << " Bcm1fSD::setDetUnitId -w- Detector name not BCM1F " << std::endl;
114 
115  int sensorNo = touch->GetReplicaNumber(0);
116  int diamondNo = touch->GetReplicaNumber(1);
117 // int detectorNo = touch->GetReplicaNumber(2);
118  int volumeNo = touch->GetReplicaNumber(3);
119 
120  // Detector ID definition
121  // detId = XYYZ
122  // X = volume, 1: +Z, 2: -Z
123  // YY = diamond, 01-12, 12: phi = 90 deg, numbering clockwise when looking from the IP
124  // Z = sensor, 1 or 2, clockwise when looking from the IP
125 
126 // detId = 10*volumeNo + diamondNo;
127  detId = 1000*volumeNo + 10*diamondNo + sensorNo;
128 
129 // for ( int ii = 0; ii < level; ii++ ) {
130 // int i = level - ii - 1;
131 // G4String name = touch->GetVolume(i)->GetName();
132 // int copyno = touch->GetReplicaNumber(i);
133 // const G4ThreeVector trans = touch->GetVolume(i)->GetTranslation();
134 // std::cout << " detId = " << detId << " name = " << name << " " << copyno
135 // << ", translation x, y, z = " << trans.x() << ", " <<trans.y() << ", " <<trans.z()
136 // << ", eta = " << trans.eta() << std::endl;
137 // }
138  }
139 
140  return detId;
141 }
142 
143 void Bcm1fSD::EndOfEvent(G4HCofThisEvent *) {
144 
145  LogDebug("Bcm1fSD")<< " Saving the last hit in a ROU " << GetName();
146 
147  if (mySimHit == nullptr) return;
148  sendHit();
149 }
150 
152  if (slave->name() == hname) { cc=slave->hits(); }
153 }
154 
156  if (mySimHit == nullptr) return;
157  LogDebug("Bcm1fSD") << " Storing PSimHit: " << pname << " " << mySimHit->detUnitId()
158  << " " << mySimHit->trackId() << " " << mySimHit->energyLoss()
159  << " " << mySimHit->entryPoint() << " " << mySimHit->exitPoint();
160 
162 
163  // clean up
164  delete mySimHit;
165  mySimHit = nullptr;
166  lastTrack = 0;
167  lastId = 0;
168 }
169 
170 void Bcm1fSD::updateHit(G4Step * aStep) {
171 
173  float theEnergyLoss = aStep->GetTotalEnergyDeposit()/GeV;
174  mySimHit->setExitPoint(theExitPoint);
175  LogDebug("Bcm1fSD")<< " Before : " << mySimHit->energyLoss();
176  mySimHit->addEnergyLoss(theEnergyLoss);
178 
179  LogDebug("Bcm1fSD") << " Updating: new exitpoint " << pname << " "
180  << theExitPoint << " new energy loss " << theEnergyLoss
181  << "\n Updated PSimHit: " << mySimHit->detUnitId()
182  << " " << mySimHit->trackId()
183  << " " << mySimHit->energyLoss() << " "
184  << mySimHit->entryPoint() << " " << mySimHit->exitPoint();
185 }
186 
187 bool Bcm1fSD::newHit(G4Step * aStep) {
188 
189  G4Track * theTrack = aStep->GetTrack();
190  uint32_t theDetUnitId = setDetUnitId(aStep);
191  unsigned int theTrackID = theTrack->GetTrackID();
192 
193  LogDebug("Bcm1fSD") << " OLD (d,t) = (" << lastId << "," << lastTrack
194  << "), new = (" << theDetUnitId << "," << theTrackID << ") return "
195  << ((theTrackID == lastTrack) && (lastId == theDetUnitId));
196  if ((mySimHit != nullptr) && (theTrackID == lastTrack) && (lastId == theDetUnitId) && closeHit(aStep))
197  return false;
198  return true;
199 }
200 
201 bool Bcm1fSD::closeHit(G4Step * aStep) {
202 
203  if (mySimHit == nullptr) return false;
204  const float tolerance = 0.05 * mm; // 50 micron are allowed between the exit
205  // point of the current hit and the entry point of the new hit
207  LogDebug("Bcm1fSD")<< " closeHit: distance = " << (mySimHit->exitPoint()-theEntryPoint).mag();
208 
209  if ((mySimHit->exitPoint()-theEntryPoint).mag()<tolerance) return true;
210  return false;
211 }
212 
213 void Bcm1fSD::createHit(G4Step * aStep) {
214 
215  if (mySimHit != nullptr) {
216  delete mySimHit;
217  mySimHit=nullptr;
218  }
219 
220  G4Track * theTrack = aStep->GetTrack();
221  G4VPhysicalVolume * v = aStep->GetPreStepPoint()->GetPhysicalVolume();
222 
225 
226  float thePabs = aStep->GetPreStepPoint()->GetMomentum().mag()/GeV;
227  float theTof = aStep->GetPreStepPoint()->GetGlobalTime()/nanosecond;
228  float theEnergyLoss = aStep->GetTotalEnergyDeposit()/GeV;
229  int theParticleType = myG4TrackToParticleID->particleID(theTrack);
230  uint32_t theDetUnitId = setDetUnitId(aStep);
231 
234  pname = theTrack->GetDynamicParticle()->GetDefinition()->GetParticleName();
235 
236  unsigned int theTrackID = theTrack->GetTrackID();
237 
238  G4ThreeVector gmd = aStep->GetPreStepPoint()->GetMomentumDirection();
239  // convert it to local frame
240  G4ThreeVector lmd = ((G4TouchableHistory *)(aStep->GetPreStepPoint()->GetTouchable()))->GetHistory()->GetTopTransform().TransformAxis(gmd);
242  float theThetaAtEntry = lnmd.theta();
243  float thePhiAtEntry = lnmd.phi();
244 
245  mySimHit = new UpdatablePSimHit(theEntryPoint,theExitPoint,thePabs,theTof,
246  theEnergyLoss,theParticleType,theDetUnitId,
247  theTrackID,theThetaAtEntry,thePhiAtEntry,
248  theG4ProcessTypeEnumerator->processId(theTrack->GetCreatorProcess()));
249 
250  LogDebug("Bcm1fSD") << " Created PSimHit: " << pname << " "
251  << mySimHit->detUnitId() << " " << mySimHit->trackId()
252  << " " << mySimHit->energyLoss() << " "
253  << mySimHit->entryPoint() << " " << mySimHit->exitPoint();
254  lastId = theDetUnitId;
255  lastTrack = theTrackID;
256  oldVolume = v;
257 }
258 
259 void Bcm1fSD::update(const BeginOfJob * ) { }
260 
262 
263  clearHits();
264  eventno = (*i)()->GetEventID();
265  mySimHit = nullptr;
266 }
267 
268 void Bcm1fSD::update(const BeginOfTrack *bot) {
269 
270  const G4Track* gTrack = (*bot)();
271  pname = gTrack->GetDynamicParticle()->GetDefinition()->GetParticleName();
272 }
273 
275  slave->Initialize();
276 }
277 
279  G4VUserTrackInformation* temp = gTrack->GetUserInformation();
280  if (temp == nullptr){
281  edm::LogError("Bcm1fSD") <<" ERROR: no G4VUserTrackInformation available";
282  abort();
283  }else{
284  TrackInformation* info = dynamic_cast<TrackInformation*>(temp);
285  if (info == nullptr){
286  edm::LogError("Bcm1fSD") <<" ERROR: TkSimTrackSelection: the UserInformation does not appear to be a TrackInformation";
287  }
288  return info;
289  }
290 }
#define LogDebug(id)
UpdatablePSimHit * mySimHit
Definition: Bcm1fSD.h:66
Local3DPoint globalExitPoint
Definition: Bcm1fSD.h:71
T getParameter(std::string const &) const
bool ProcessHits(G4Step *, G4TouchableHistory *) override
Definition: Bcm1fSD.cc:60
static const TGPicture * info(bool iBackgroundIsBlack)
const double GeV
Definition: MathUtil.h:16
G4TrackToParticleID * myG4TrackToParticleID
Definition: Bcm1fSD.h:65
uint32_t lastId
Definition: Bcm1fSD.h:73
bool storeTrack() const
Local3DPoint ConvertToLocal3DPoint(const G4ThreeVector &point) const
G4VPhysicalVolume * oldVolume
Definition: Bcm1fSD.h:72
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
void clearHits() override
Definition: Bcm1fSD.cc:274
std::string name() const
G4ProcessTypeEnumerator * theG4ProcessTypeEnumerator
Definition: Bcm1fSD.h:64
TrackInformation * getOrCreateTrackInformation(const G4Track *)
Definition: Bcm1fSD.cc:278
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
void setExitPoint(const Local3DPoint &exit)
~Bcm1fSD() override
Definition: Bcm1fSD.cc:54
virtual bool closeHit(G4Step *)
Definition: Bcm1fSD.cc:201
TrackingSlaveSD * slave
Definition: Bcm1fSD.h:63
Local3DPoint globalEntryPoint
Definition: Bcm1fSD.h:70
#define nullptr
float energyHistoryCut
Definition: Bcm1fSD.h:68
type of data representation of DDCompactView
Definition: DDCompactView.h:90
int eventno
Definition: Bcm1fSD.h:75
virtual bool newHit(G4Step *)
Definition: Bcm1fSD.cc:187
uint32_t setDetUnitId(const G4Step *) override
Definition: Bcm1fSD.cc:91
Geom::Theta< T > theta() const
Definition: PV3DBase.h:75
std::vector< PSimHit > & hits()
unsigned int processId(const G4VProcess *p) const
Local3DPoint exitPoint() const
Exit point in the local Det frame.
Definition: PSimHit.h:38
void EndOfEvent(G4HCofThisEvent *) override
Definition: Bcm1fSD.cc:143
virtual void Initialize()
Local3DPoint FinalStepPosition(const G4Step *step, coordinates) const
float energyCut
Definition: Bcm1fSD.h:67
static int particleID(const G4Track *)
void addEnergyLoss(float eloss)
Bcm1fSD(const std::string &, const DDCompactView &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
Definition: Bcm1fSD.cc:33
void fillHits(edm::PSimHitContainer &, const std::string &) override
Definition: Bcm1fSD.cc:151
float energyLoss() const
The energy deposit in the PSimHit, in ???.
Definition: PSimHit.h:75
virtual void sendHit()
Definition: Bcm1fSD.cc:155
virtual void updateHit(G4Step *)
Definition: Bcm1fSD.cc:170
void update(const BeginOfEvent *) override
This routine will be called when the appropriate signal arrives.
Definition: Bcm1fSD.cc:261
unsigned int trackId() const
Definition: PSimHit.h:102
unsigned int lastTrack
Definition: Bcm1fSD.h:74
virtual bool processHits(const PSimHit &)
std::vector< PSimHit > PSimHitContainer
virtual void createHit(G4Step *)
Definition: Bcm1fSD.cc:213
std::string pname
Definition: Bcm1fSD.h:76
Local3DPoint entryPoint() const
Entry point in the local Det frame.
Definition: PSimHit.h:35
unsigned int detUnitId() const
Definition: PSimHit.h:93
Local3DPoint InitialStepPosition(const G4Step *step, coordinates) const