CMS 3D CMS Logo

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