CMS 3D CMS Logo

TotemSD.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: Forward
4 // Class : TotemSD
5 //
6 // Implementation:
7 // <Notes on implementation>
8 //
9 // Original Author:
10 // Created: Tue May 16 10:14:34 CEST 2006
11 //
12 
13 // system include files
14 
15 // user include files
19 
23 
26 
32 
35 
36 #include "G4SDManager.hh"
37 #include "G4Step.hh"
38 #include "G4Track.hh"
39 #include "G4VProcess.hh"
40 
41 #include "G4PhysicalConstants.hh"
42 #include "G4SystemOfUnits.hh"
43 
44 
45 //
46 // constructors and destructor
47 //
49  const SensitiveDetectorCatalog & clg,
50  edm::ParameterSet const & p, const SimTrackManager* manager) :
51  SensitiveTkDetector(name, cpv, clg, p), numberingScheme(nullptr),
52  hcID(-1), theHC(nullptr), theManager(manager), currentHit(nullptr), theTrack(nullptr),
53  currentPV(nullptr), unitID(0), previousUnitID(0), preStepPoint(nullptr),
54  postStepPoint(nullptr){
55 
56  //Parameters
58  int verbn = m_p.getUntrackedParameter<int>("Verbosity");
59 
60  SetVerboseLevel(verbn);
61 
62  slave = new TrackingSlaveSD(name);
63 
64  if (name == "TotemHitsT1") {
66  } else if (name == "TotemHitsT2Si") {
68  } else if (name == "TotemHitsT2Gem") {
70  } else if (name == "TotemHitsRP") {
72  } else {
73  edm::LogWarning("ForwardSim") << "TotemSD: ReadoutName not supported\n";
74  }
75 }
76 
78  delete slave;
79  delete numberingScheme;
80 }
81 
82 bool TotemSD::ProcessHits(G4Step * aStep, G4TouchableHistory * ) {
83 
84  getStepInfo(aStep);
85  if (!hitExists() && edeposit>0.) {
86  createNewHit();
87  return true;
88  }
89  if (!hitExists() && (((unitID==1111 || unitID==2222) &&
90  ParentId==0 && ParticleType==2212))) {
92  }
93  return true;
94 }
95 
96 uint32_t TotemSD::setDetUnitId(const G4Step * aStep) {
97 
98  return (numberingScheme == nullptr ? 0 : numberingScheme->getUnitID(aStep));
99 }
100 
101 void TotemSD::Initialize(G4HCofThisEvent * HCE) {
102 
103  LogDebug("ForwardSim") << "TotemSD : Initialize called for " << GetName();
104 
105  theHC = new TotemG4HitCollection(GetName(), collectionName[0]);
106  if (hcID<0)
107  hcID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
108  HCE->AddHitsCollection(hcID, theHC);
109 
110  tsID = -2;
111  primID = -2;
112 }
113 
114 void TotemSD::EndOfEvent(G4HCofThisEvent* ) {
115 
116  // here we loop over transient hits and make them persistent
117  for (int j=0; j<theHC->entries() && j<15000; j++) {
118  TotemG4Hit* aHit = (*theHC)[j];
119 #ifdef ddebug
120  LogDebug("ForwardSim") << "HIT NUMERO " << j << "unit ID = "
121  << aHit->getUnitID() << "\n"
122  << " " << "enrty z "
123  << aHit->getEntry().z() << "\n"
124  << " " << "theta "
125  << aHit->getThetaAtEntry() << "\n";
126 #endif
127  Local3DPoint theExitPoint(0,0,0);
128  Local3DPoint Entrata(aHit->getEntry().x(),
129  aHit->getEntry().y(),
130  aHit->getEntry().z());
131  slave->processHits(PSimHit(Entrata,theExitPoint,
132  aHit->getPabs(), aHit->getTof(),
133  aHit->getEnergyLoss(), aHit->getParticleType(),
134  aHit->getUnitID(), aHit->getTrackID(),
135  aHit->getThetaAtEntry(),aHit->getPhiAtEntry()));
136 
137  }
138 }
139 
141  LogDebug("ForwardSim") << "TotemSD: Collection " << theHC->GetName();
142  theHC->PrintAllHits();
143 }
144 
146  if (slave->name() == hname) { cc=slave->hits(); }
147 }
148 
150  LogDebug("ForwardSim") << " Dispatched BeginOfEvent for " << GetName()
151  << " !" ;
152  clearHits();
153 }
154 
156  slave->Initialize();
157 }
158 
159 G4ThreeVector TotemSD::setToLocal(const G4ThreeVector& global) {
160 
161  G4ThreeVector localPoint;
162  const G4VTouchable* touch= preStepPoint->GetTouchable();
163  localPoint = touch->GetHistory()->GetTopTransform().TransformPoint(global);
164  return localPoint;
165 }
166 
167 void TotemSD::getStepInfo(const G4Step* aStep) {
168 
169  preStepPoint = aStep->GetPreStepPoint();
170  postStepPoint= aStep->GetPostStepPoint();
171  theTrack = aStep->GetTrack();
172  hitPoint = preStepPoint->GetPosition();
173  currentPV = preStepPoint->GetPhysicalVolume();
174 
175  // double weight = 1;
176  G4String name = currentPV->GetName();
177  name.assign(name,0,4);
178  G4String particleType = theTrack->GetDefinition()->GetParticleName();
179  edeposit = aStep->GetTotalEnergyDeposit();
180 
181  tSlice = (postStepPoint->GetGlobalTime() )/nanosecond;
182  tSliceID = (int) tSlice;
183  unitID = setDetUnitId(aStep);
184 #ifdef debug
185  LogDebug("ForwardSim") << "UNITa " << unitID;
186 #endif
187  primaryID = theTrack->GetTrackID();
188 
189  Posizio = hitPoint;
190  Pabs = aStep->GetPreStepPoint()->GetMomentum().mag()/GeV;
191  Tof = aStep->GetPostStepPoint()->GetGlobalTime()/nanosecond;
192 
193  Eloss = aStep->GetTotalEnergyDeposit()/GeV;
194  ParticleType = theTrack->GetDefinition()->GetPDGEncoding();
195 
196  ThetaAtEntry = aStep->GetPreStepPoint()->GetPosition().theta()/deg;
197  PhiAtEntry = aStep->GetPreStepPoint()->GetPosition().phi()/deg;
198 
199  ParentId = theTrack->GetParentID();
200  Vx = theTrack->GetVertexPosition().x();
201  Vy = theTrack->GetVertexPosition().y();
202  Vz = theTrack->GetVertexPosition().z();
203 }
204 
206 
207  if (primaryID<1) {
208  edm::LogWarning("ForwardSim") << "***** TotemSD error: primaryID = "
209  << primaryID
210  << " maybe detector name changed";
211  }
212 
213  // Update if in the same detector, time-slice and for same track
214  // if (primaryID == primID && tSliceID == tsID && unitID==previousUnitID) {
215  if (tSliceID == tsID && unitID==previousUnitID) {
216  updateHit();
217  return true;
218  }
219 
220  // Reset entry point for new primary
221  if (primaryID != primID)
223 
224  //look in the HitContainer whether a hit with the same primID, unitID,
225  //tSliceID already exists:
226 
227  bool found = false;
228 
229  for (int j=0; j<theHC->entries()&&!found; j++) {
230  TotemG4Hit* aPreviousHit = (*theHC)[j];
231  if (aPreviousHit->getTrackID() == primaryID &&
232  aPreviousHit->getTimeSliceID() == tSliceID &&
233  aPreviousHit->getUnitID() == unitID ) {
234  currentHit = aPreviousHit;
235  found = true;
236  break;
237  }
238  }
239 
240  if (found) {
241  updateHit();
242  return true;
243  } else {
244  return false;
245  }
246 }
247 
249 
250 #ifdef debug
251  LogDebug("ForwardSim") << "TotemSD CreateNewHit for"
252  << " PV " << currentPV->GetName()
253  << " PVid = " << currentPV->GetCopyNo()
254  << " MVid = " << currentPV->GetMother()->GetCopyNo()
255  << " Unit " << unitID << "\n"
256  << " primary " << primaryID
257  << " time slice " << tSliceID
258  << " For Track " << theTrack->GetTrackID()
259  << " which is a "
260  << theTrack->GetDefinition()->GetParticleName();
261 
262  if (theTrack->GetTrackID()==1) {
263  LogDebug("ForwardSim") << " of energy " << theTrack->GetTotalEnergy();
264  } else {
265  LogDebug("ForwardSim") << " daughter of part. " << theTrack->GetParentID();
266  }
267 
268  if (theTrack->GetCreatorProcess()!=nullptr)
269  LogDebug("ForwardSim") << theTrack->GetCreatorProcess()->GetProcessName() ;
270  else
271  LogDebug("ForwardSim") << "NO process";
272 #endif
273 
274 
275  currentHit = new TotemG4Hit;
280 
287 
288  currentHit->setEntry(Posizio.x(),Posizio.y(),Posizio.z());
289 
291  currentHit->setVx(Vx);
292  currentHit->setVy(Vy);
293  currentHit->setVz(Vz);
294 
295  updateHit();
296 
298 }
299 
301 
302 // LogDebug("ForwardSim") << "INSIDE CREATE NEW HIT EVO ";
303 
304  currentHit = new TotemG4Hit;
309 
316 
317  // LogDebug("ForwardSim") << Posizio.x() << " " << Posizio.y() << " " << Posizio.z();
318 
320  currentHit->setVx(Vx);
321  currentHit->setVy(Vy);
322  currentHit->setVz(Vz);
323 
324  G4ThreeVector _PosizioEvo;
325  int flagAcc=0;
326  _PosizioEvo=posizioEvo(Posizio,Vx,Vy,Vz,Pabs,flagAcc);
327 
328  if(flagAcc==1){
329  currentHit->setEntry(_PosizioEvo.x(),_PosizioEvo.y(),_PosizioEvo.z());
330 
331  updateHit();
332 
334  }
335  // LogDebug("ForwardSim") << "STORED HIT IN: " << unitID;
336 }
337 
338 G4ThreeVector TotemSD::posizioEvo(const G4ThreeVector& Pos, double vx, double vy,
339  double vz, double pabs, int& accettanza) {
340  accettanza=0;
341  //Pos.xyz() in mm
342  G4ThreeVector PosEvo;
343  double ThetaX=std::atan((Pos.x()-vx)/(Pos.z()-vz));
344  double ThetaY=std::atan((Pos.y()-vy)/(Pos.z()-vz));
345  double X_at_0 =(vx-((Pos.x()-vx)/(Pos.z()-vz))*vz)/1000.;
346  double Y_at_0 =(vy-((Pos.y()-vy)/(Pos.z()-vz))*vz)/1000.;
347 
348  // double temp_evo_X;
349  // double temp_evo_Y;
350  // double temp_evo_Z;
351  // temp_evo_Z = std::abs(Pos.z())/Pos.z()*220000.;
352 
353  //csi=-dp/d
354  double csi = std::abs((7000.-pabs)/7000.);
355 
356  // all in m
357  const int no_rp=4;
358  double x_par[no_rp+1];
359  double y_par[no_rp+1];
360  //rp z position
361  double rp[no_rp]={141.,149.,198.,220.};
362  //{lx0,mlx} for each rp; Lx=lx0+mlx*csi
363  double leffx[][2]={{122.5429,-46.9312},{125.4194,-49.1849},{152.6,-81.157},{98.8914,-131.8390}};
364  //{ly0,mly} for each rp; Ly=ly0+mly*csi
365  double leffy[][2]={{124.2314,-55.4852},{127.7825,-57.4503},{179.455,-76.274},{273.0931,-40.4626}};
366  //{vx0,mvx0} for each rp; vx=vx0+mvx*csi
367  double avx[][2]={{0.515483,-1.0123},{0.494122,-1.0534},{0.2217,-1.483},{0.004633,-1.0719}};
368  //{vy0,mvy0} for each rp; vy=vy0+mvy*csi
369  double avy[][2]={{0.371418,-1.6327},{0.349035,-1.6955},{0.0815,-2.59},{0.007592,-4.0841}};
370  //{D0,md,a,b} for each rp; D=D0+(md+a*thetax)*csi+b*thetax
371  double ddx[][4]= {{-0.082336,-0.092513,112.3436,-82.5029},{-0.086927,-0.097670,114.9513,-82.9835},
372  {-0.092117,-0.0915,180.6236,-82.443},{-0.050470,0.058837,208.1106,20.8198}};
373  // {10sigma_x+0.5mm,10sigma_y+0.5mm}
374  double detlim[][2]={{0,0},{0.0028,0.0021},{0,0},{0.0008,0.0013}};
375  //{rmax,dmax}
376  double pipelim[][2]={{0.026,0.026},{0.04,0.04},{0.0226,0.0177},{0.04,0.04}};
377 
378 
379  for(int j=0; j<no_rp ; j++) {
380  //y=Ly*thetay+vy*y0
381  //x=Lx*thetax+vx*x0-csi*D
382  y_par[j]=ThetaY*(leffy[j][0]+leffy[j][1]*csi)+(avy[j][0]+avy[j][1]*csi)*Y_at_0;
383  x_par[j]=ThetaX*(leffx[j][0]+leffx[j][1]*csi)+(avx[j][0]+avx[j][1]*csi)*X_at_0-
384  csi*(ddx[j][0]+(ddx[j][1]+ddx[j][2]*ThetaX)*csi+ddx[j][3]*ThetaX);
385  }
386 
387 
388  //pass TAN@141
389  if (std::abs(y_par[0])<pipelim[0][1] && sqrt((y_par[0]*y_par[0])+(x_par[0]*x_par[0]))<pipelim[0][0]) {
390  //pass 149
391  if ((sqrt((y_par[1]*y_par[1])+(x_par[1]*x_par[1]))<pipelim[1][0]) &&
392  (std::abs(y_par[1])>detlim[1][1] || x_par[1]>detlim[1][0])) {
393  accettanza = 1;
394  }
395  }
396 
397 
398  //pass TAN@141
399  if (std::abs(y_par[0])<pipelim[0][1] && sqrt((y_par[0])*(y_par[0])+(x_par[0])*(x_par[0]))<pipelim[0][0]) {
400  //pass Q5@198
401  if (std::abs(y_par[2])<pipelim[2][1] && sqrt((y_par[2]*y_par[2])+(x_par[2]*x_par[2]))<pipelim[2][0]) {
402  //pass 220
403  if ((sqrt((y_par[3]*y_par[3])+(x_par[3]*x_par[3]))<pipelim[3][0]) &&
404  (std::abs(y_par[3])>detlim[3][1] || x_par[3]>detlim[3][0])) {
405  accettanza = 1;
406 
407  PosEvo.setX(1000*x_par[3]);
408  PosEvo.setY(1000*y_par[3]);
409  PosEvo.setZ(1000*rp[3]);
410  if(Pos.z()<vz)PosEvo.setZ(-1000*rp[3]);
411  }
412  }
413 
414  }
415 /*
416  LogDebug("ForwardSim") << "\n"
417  << "ACCETTANZA: "<<accettanza << "\n"
418  << "CSI: "<< csi << "\n"
419  << "Theta_X: " << ThetaX << "\n"
420  << "Theta_Y: " << ThetaY << "\n"
421  << "X_at_0: "<< X_at_0 << "\n"
422  << "Y_at_0: "<< Y_at_0 << "\n"
423  << "x_par[3]: "<< x_par[3] << "\n"
424  << "y_par[3]: "<< y_par[3] << "\n"
425  << "pos " << Pos.x() << " " << Pos.y() << " "
426  << Pos.z() << "\n" << "V "<< vx << " " << vy << " "
427  << vz << "\n"
428 */
429 // --------------
430  return PosEvo;
431 }
432 
434  //
435  if (Eloss > 0.) {
436 
437 #ifdef debug
438  LogDebug("ForwardSim") << "G4TotemT1SD updateHit: add eloss " << Eloss
439  << "\nCurrentHit=" << currentHit
440  << ", PostStepPoint="
441  << postStepPoint->GetPosition();
442 #endif
443 
445  }
446 
447  // buffer for next steps:
448  tsID = tSliceID;
449  primID = primaryID;
451 }
452 
454 
455  if (primID<0) return;
456  if (hit == nullptr ) {
457  edm::LogWarning("ForwardSim") << "TotemSD: hit to be stored is NULL !!";
458  return;
459  }
460 
461  theHC->insert( hit );
462 }
463 
465 
467  incidentEnergy = preStepPoint->GetKineticEnergy();
468 }
469 
470 
#define LogDebug(id)
void setVz(float p)
Definition: TotemG4Hit.cc:183
G4ThreeVector setToLocal(const G4ThreeVector &globalPoint)
Definition: TotemSD.cc:159
void resetForNewPrimary()
Definition: TotemSD.cc:464
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
int ParentId
Definition: TotemSD.h:118
void Initialize(G4HCofThisEvent *HCE) override
Definition: TotemSD.cc:101
void setPabs(float e)
Definition: TotemG4Hit.cc:153
uint32_t unitID
Definition: TotemSD.h:100
const double GeV
Definition: MathUtil.h:16
float Pabs
Definition: TotemSD.h:110
G4ThreeVector posizioEvo(const G4ThreeVector &, double, double, double, double, int &)
Definition: TotemSD.cc:338
G4ThreeVector Posizio
Definition: TotemSD.h:109
TotemG4HitCollection * theHC
Definition: TotemSD.h:93
float Tof
Definition: TotemSD.h:111
void setEnergyLoss(float e)
Definition: TotemG4Hit.cc:155
G4ThreeVector entrancePoint
Definition: TotemSD.h:88
void setEntry(double x, double y, double z)
Definition: TotemG4Hit.h:57
double tSlice
Definition: TotemSD.h:102
float getTof() const
Definition: TotemG4Hit.cc:149
int primaryID
Definition: TotemSD.h:101
bool ProcessHits(G4Step *, G4TouchableHistory *) override
Definition: TotemSD.cc:82
void createNewHit()
Definition: TotemSD.cc:248
std::string name() const
uint32_t previousUnitID
Definition: TotemSD.h:100
float getEnergyLoss() const
Definition: TotemG4Hit.cc:150
float incidentEnergy
Definition: TotemSD.h:89
void storeHit(TotemG4Hit *)
Definition: TotemSD.cc:453
#define nullptr
void update(const BeginOfEvent *) override
This routine will be called when the appropriate signal arrives.
Definition: TotemSD.cc:149
math::XYZPoint getEntry() const
Definition: TotemG4Hit.cc:123
TotemVDetectorOrganization * numberingScheme
Definition: TotemSD.h:81
G4int hcID
Definition: TotemSD.h:92
int tSliceID
Definition: TotemSD.h:101
float getPabs() const
Definition: TotemG4Hit.cc:148
Compact representation of the geometrical detector hierarchy.
Definition: DDCompactView.h:80
void setVy(float p)
Definition: TotemG4Hit.cc:180
static TrackerG4SimHitNumberingScheme & numberingScheme(const DDCompactView &cpv, const GeometricDet &det)
const G4StepPoint * preStepPoint
Definition: TotemSD.h:104
void setTimeSlice(double d)
Definition: TotemG4Hit.cc:141
std::vector< PSimHit > & hits()
std::string const collectionName[nCollections]
Definition: Collections.h:47
virtual void Initialize()
float Eloss
Definition: TotemSD.h:112
float getPhiAtEntry() const
Definition: TotemG4Hit.cc:159
T sqrt(T t)
Definition: SSEVec.h:18
float edeposit
Definition: TotemSD.h:106
float Vx
Definition: TotemSD.h:119
G4ThreeVector hitPoint
Definition: TotemSD.h:107
void setVx(float p)
Definition: TotemG4Hit.cc:177
void fillHits(edm::PSimHitContainer &, const std::string &) override
Definition: TotemSD.cc:145
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int getTrackID() const
Definition: TotemG4Hit.cc:134
int getTimeSliceID() const
Definition: TotemG4Hit.cc:142
~TotemSD() override
Definition: TotemSD.cc:77
void setParentId(int p)
Definition: TotemG4Hit.cc:174
TotemG4Hit * currentHit
Definition: TotemSD.h:97
G4VPhysicalVolume * currentPV
Definition: TotemSD.h:99
void setIncidentEnergy(double e)
Definition: TotemG4Hit.cc:132
uint32_t setDetUnitId(const G4Step *) override
Definition: TotemSD.cc:96
void PrintAll() override
Definition: TotemSD.cc:140
short ParticleType
Definition: TotemSD.h:113
void setThetaAtEntry(float t)
Definition: TotemG4Hit.cc:161
float PhiAtEntry
Definition: TotemSD.h:116
TotemSD(const std::string &, const DDCompactView &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
Definition: TotemSD.cc:48
const G4StepPoint * postStepPoint
Definition: TotemSD.h:105
void setUnitID(uint32_t i)
Definition: TotemG4Hit.cc:138
void setTrackID(int i)
Definition: TotemG4Hit.cc:135
void setPhiAtEntry(float f)
Definition: TotemG4Hit.cc:162
void EndOfEvent(G4HCofThisEvent *eventHC) override
Definition: TotemSD.cc:114
uint32_t getUnitID() const
Definition: TotemG4Hit.cc:137
void updateHit()
Definition: TotemSD.cc:433
virtual uint32_t getUnitID(const G4Step *aStep) const =0
void setParticleType(short i)
Definition: TotemG4Hit.cc:156
void setTof(float e)
Definition: TotemG4Hit.cc:154
float Vz
Definition: TotemSD.h:119
G4int primID
Definition: TotemSD.h:90
TrackingSlaveSD * slave
Definition: TotemSD.h:80
void createNewHitEvo()
Definition: TotemSD.cc:300
int getParticleType() const
Definition: TotemG4Hit.cc:151
G4Track * theTrack
Definition: TotemSD.h:98
virtual bool processHits(const PSimHit &)
std::vector< PSimHit > PSimHitContainer
float getThetaAtEntry() const
Definition: TotemG4Hit.cc:158
int tsID
Definition: TotemSD.h:96
G4THitsCollection< TotemG4Hit > TotemG4HitCollection
float Vy
Definition: TotemSD.h:119
void clearHits() override
Definition: TotemSD.cc:155
float ThetaAtEntry
Definition: TotemSD.h:115
bool hitExists()
Definition: TotemSD.cc:205
void getStepInfo(const G4Step *aStep)
Definition: TotemSD.cc:167