CMS 3D CMS Logo

CastorSD.cc
Go to the documentation of this file.
1 // File: CastorSD.cc
3 // Date: 02.04
4 // UpDate: 07.04 - C3TF & C4TF semi-trapezoid added
5 // Description: Sensitive Detector class for Castor
7 
10 
13 
14 #include "G4SDManager.hh"
15 #include "G4Step.hh"
16 #include "G4Track.hh"
17 #include "G4VProcess.hh"
18 
19 #include "G4ios.hh"
20 #include "G4Cerenkov.hh"
21 #include "G4LogicalVolumeStore.hh"
22 
23 #include "CLHEP/Units/GlobalSystemOfUnits.h"
24 #include "CLHEP/Units/GlobalPhysicalConstants.h"
25 #include "Randomize.hh"
26 #include "G4Poisson.hh"
27 
28 //#define debugLog
29 
31  const SensitiveDetectorCatalog & clg,
32  edm::ParameterSet const & p,
33  const SimTrackManager* manager) :
34  CaloSD(name, cpv, clg, p, manager), numberingScheme(nullptr), lvC3EF(nullptr),
35  lvC3HF(nullptr), lvC4EF(nullptr), lvC4HF(nullptr), lvCAST(nullptr) {
36 
37  edm::ParameterSet m_CastorSD = p.getParameter<edm::ParameterSet>("CastorSD");
38  useShowerLibrary = m_CastorSD.getParameter<bool>("useShowerLibrary");
39  energyThresholdSL = m_CastorSD.getParameter<double>("minEnergyInGeVforUsingSLibrary");
40  energyThresholdSL = energyThresholdSL*GeV; // Convert GeV => MeV
41 
42  non_compensation_factor = m_CastorSD.getParameter<double>("nonCompensationFactor");
43 
44  if (useShowerLibrary) showerLibrary = new CastorShowerLibrary(name, p);
45 
47 
48  edm::LogInfo("ForwardSim")
49  << "***************************************************\n"
50  << "* *\n"
51  << "* Constructing a CastorSD with name " << GetName() << "\n"
52  << "* *\n"
53  << "***************************************************";
54 
55  const G4LogicalVolumeStore * lvs = G4LogicalVolumeStore::GetInstance();
56  std::vector<G4LogicalVolume*>::const_iterator lvcite;
57  for (lvcite = lvs->begin(); lvcite != lvs->end(); lvcite++) {
58  if (strcmp(((*lvcite)->GetName()).c_str(),"C3EF") == 0) lvC3EF = (*lvcite);
59  if (strcmp(((*lvcite)->GetName()).c_str(),"C3HF") == 0) lvC3HF = (*lvcite);
60  if (strcmp(((*lvcite)->GetName()).c_str(),"C4EF") == 0) lvC4EF = (*lvcite);
61  if (strcmp(((*lvcite)->GetName()).c_str(),"C4HF") == 0) lvC4HF = (*lvcite);
62  if (strcmp(((*lvcite)->GetName()).c_str(),"CAST") == 0) lvCAST = (*lvcite);
63  if (lvC3EF != nullptr && lvC3HF != nullptr && lvC4EF != nullptr && lvC4HF != nullptr && lvCAST != nullptr) break;
64  }
65  edm::LogInfo("ForwardSim") << "CastorSD:: LogicalVolume pointers\n"
66  << lvC3EF << " for C3EF; " << lvC3HF
67  << " for C3HF; " << lvC4EF << " for C4EF; "
68  << lvC4HF << " for C4HF; "
69  << lvCAST << " for CAST. " << std::endl;
70 
71  // if(useShowerLibrary) edm::LogInfo("ForwardSim") << "\n Using Castor Shower Library \n";
72 
73 }
74 
75 //=============================================================================================
76 
78  if (useShowerLibrary) delete showerLibrary;
79 }
80 
81 //=============================================================================================
82 
84  if (useShowerLibrary) {
85  // showerLibrary = new CastorShowerLibrary(name, cpv, p);
86  G4ParticleTable * theParticleTable = G4ParticleTable::GetParticleTable();
87  showerLibrary->initParticleTable(theParticleTable);
88  edm::LogInfo("ForwardSim") << "CastorSD::initRun: Using Castor Shower Library \n";
89  }
90 }
91 
92 //=============================================================================================
93 
94 double CastorSD::getEnergyDeposit(G4Step * aStep) {
95 
96  double NCherPhot = 0.;
97 
98  // Get theTrack
99  G4Track* theTrack = aStep->GetTrack();
100 
101  // preStepPoint information *********************************************
102 
103  G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
104  G4VPhysicalVolume* currentPV = preStepPoint->GetPhysicalVolume();
105  G4LogicalVolume* currentLV = currentPV->GetLogicalVolume();
106 
107 #ifdef debugLog
108  G4String name = currentPV->GetName();
109  std::string nameVolume;
110  nameVolume.assign(name,0,4);
111 
112  G4SteppingControl stepControlFlag = aStep->GetControlFlag();
113  if (aStep->IsFirstStepInVolume())
114  LogDebug("ForwardSim") << "CastorSD::getEnergyDeposit:"
115  << "\n IsFirstStepInVolume " ;
116 #endif
117 
118 
119 
120 #ifdef debugLog
121  if (useShowerLibrary && currentLV==lvCAST) {
122  LogDebug("ForwardSim") << "CastorSD::getEnergyDeposit:"
123  << "\n TrackID , ParentID , ParticleName ,"
124  << " eta , phi , z , time ,"
125  << " K , E , Mom "
126  << "\n TRACKINFO: "
127  << theTrack->GetTrackID()
128  << " , "
129  << theTrack->GetParentID()
130  << " , "
131  << theTrack->GetDefinition()->GetParticleName()
132  << " , "
133  << theTrack->GetPosition().eta()
134  << " , "
135  << theTrack->GetPosition().phi()
136  << " , "
137  << theTrack->GetPosition().z()
138  << " , "
139  << theTrack->GetGlobalTime()
140  << " , "
141  << theTrack->GetKineticEnergy()
142  << " , "
143  << theTrack->GetTotalEnergy()
144  << " , "
145  << theTrack->GetMomentum().mag() ;
146  if(theTrack->GetTrackID() != 1)
147  LogDebug("ForwardSim") << "CastorSD::getEnergyDeposit:"
148  << "\n CurrentStepNumber , TrackID , Particle , VertexPosition ,"
149  << " LogicalVolumeAtVertex , CreatorProcess"
150  << "\n TRACKINFO2: "
151  << theTrack->GetCurrentStepNumber()
152  << " , "
153  << theTrack->GetTrackID()
154  << " , "
155  << theTrack->GetDefinition()->GetParticleName()
156  << " , "
157  << theTrack->GetVertexPosition()
158  << " , "
159  << theTrack->GetLogicalVolumeAtVertex()->GetName()
160  << " , "
161  << theTrack->GetCreatorProcess()->GetProcessName() ;
162  } // end of if(useShowerLibrary)
163 #endif
164 
165  // if particle moves from interaction point or "backwards (halo)
166  bool backward = false;
167  const G4ThreeVector& hitPoint = preStepPoint->GetPosition();
168  const G4ThreeVector& hit_mom = preStepPoint->GetMomentumDirection();
169  double zint = hitPoint.z();
170  double pz = hit_mom.z();
171 
172  // Check if theTrack moves backward
173  if (pz * zint < 0.) backward = true;
174 
175  // Check that theTrack is above the energy threshold to use Shower Library
176  bool aboveThreshold = false;
177  if(theTrack->GetKineticEnergy() > energyThresholdSL) aboveThreshold = true;
178 
179  // Check if theTrack is a muon (if so, DO NOT use Shower Library)
180  bool notaMuon = true;
181  G4int mumPDG = 13;
182  G4int mupPDG = -13;
183  G4int parCode = theTrack->GetDefinition()->GetPDGEncoding();
184  if (parCode == mupPDG || parCode == mumPDG ) notaMuon = false;
185 
186  // angle condition
187  double theta_max = M_PI - 3.1305; // angle in radians corresponding to -5.2 eta
188  double R_mom=sqrt(hit_mom.x()*hit_mom.x() + hit_mom.y()*hit_mom.y());
189  double theta = atan2(R_mom,std::abs(pz));
190  bool angleok = false;
191  if ( theta < theta_max) angleok = true;
192 
193  // OkToUse
194  double R = sqrt(hitPoint.x()*hitPoint.x() + hitPoint.y()*hitPoint.y());
195  bool dot = false;
196  if ( zint < -14450. && R < 45.) dot = true;
197  bool inRange = true;
198  if ( zint < -14700. || R > 193.) inRange = false;
199  bool OkToUse = false;
200  if ( inRange && !dot) OkToUse = true;
201 
202  const bool particleWithinShowerLibrary = aboveThreshold &&
203  notaMuon && (!backward) && OkToUse && angleok && currentLV == lvCAST;
204 
205  if (useShowerLibrary && particleWithinShowerLibrary) {
206  // Use Castor shower library if energy is above threshold, is not a muon
207  // and is not moving backward
208  getFromLibrary(aStep);
209 
210 #ifdef debugLog
211  LogDebug("ForwardSim") << " Current logical volume is " << nameVolume ;
212 #endif
213 
214  // track is killed in getFromLibrary...
215 
216  return 0;
217  }
218 
219 
220  // Full - Simulation starts here
221 
222  double meanNCherPhot=0;
223 
224  // remember primary particle hitting the CASTOR detector
225 
226  TrackInformationExtractor TIextractor;
227  TrackInformation& trkInfo = TIextractor(theTrack);
228  if (!trkInfo.hasCastorHit()) {
229  trkInfo.setCastorHitPID(parCode);
230  }
231  const int castorHitPID = trkInfo.getCastorHitPID();
232 
233  // Check whether castor hit track is HAD
234  const bool isHad = !(castorHitPID==emPDG || castorHitPID==epPDG || castorHitPID==gammaPDG || castorHitPID == mupPDG || castorHitPID == mumPDG);
235 
236 
237  // Usual calculations
238  // G4ThreeVector hitPoint = preStepPoint->GetPosition();
239  // G4ThreeVector hit_mom = preStepPoint->GetMomentumDirection();
240  G4double stepl = aStep->GetStepLength()/cm;
241  G4double beta = preStepPoint->GetBeta();
242  G4double charge = preStepPoint->GetCharge();
243  // G4VProcess* curprocess = preStepPoint->GetProcessDefinedStep();
244  // G4String namePr = preStepPoint->GetProcessDefinedStep()->GetProcessName();
245  // std::string nameProcess;
246  // nameProcess.assign(namePr,0,4);
247 
248  // G4LogicalVolume* lv = currentPV->GetLogicalVolume();
249  // G4Material* mat = lv->GetMaterial();
250  // G4double rad = mat->GetRadlen();
251 
252 
253 #ifdef debugLog
254  // postStepPoint information *********************************************
255  G4StepPoint* postStepPoint= aStep->GetPostStepPoint();
256  G4VPhysicalVolume* postPV= postStepPoint->GetPhysicalVolume();
257 
258  G4String postname = postPV->GetName();
259  std::string postnameVolume;
260  postnameVolume.assign(postname,0,4);
261 
262  // theTrack information *************************************************
263  // G4Track* theTrack = aStep->GetTrack();
264  //G4double entot = theTrack->GetTotalEnergy();
265  G4ThreeVector vert_mom = theTrack->GetVertexMomentumDirection();
266 
267  G4ThreeVector localPoint = theTrack->GetTouchable()->GetHistory()->
268  GetTopTransform().TransformPoint(hitPoint);
269 
270  G4String particleType = theTrack->GetDefinition()->GetParticleName();
271 
272  // calculations... *************************************************
273  double phi = -100.;
274  if (vert_mom.x() != 0) phi = atan2(vert_mom.y(),vert_mom.x());
275  if (phi < 0.) phi += twopi;
276 
277 
278  double costheta =vert_mom.z()/sqrt(vert_mom.x()*vert_mom.x()+
279  vert_mom.y()*vert_mom.y()+
280  vert_mom.z()*vert_mom.z());
281  double theta = acos(std::min(std::max(costheta,double(-1.)),double(1.)));
282  double eta = -log(tan(theta/2));
283  G4int primaryID = theTrack->GetTrackID();
284  // *************************************************
285 
286 
287  // *************************************************
288  double edep = aStep->GetTotalEnergyDeposit();
289 #endif
290 
291  // *************************************************
292 
293 
294  // *************************************************
295  // take into account light collection curve for plate
296  // double weight = curve_Castor(nameVolume, preStepPoint);
297  // double edep = aStep->GetTotalEnergyDeposit() * weight;
298  // *************************************************
299 
300 
301  // *************************************************
302  /* comments for sensitive volumes:
303  C001 ...-... CP06,CBU1 ...-...CALM --- > fibres and bundle
304  for first release of CASTOR
305  CASF --- > quartz plate for first and second releases of CASTOR
306  GF2Q, GFNQ, GR2Q, GRNQ
307  for tests with my own test geometry of HF (on ask of Gavrilov)
308  C3TF, C4TF - for third release of CASTOR
309  */
310  if (currentLV == lvC3EF || currentLV == lvC4EF || currentLV == lvC3HF ||
311  currentLV == lvC4HF) {
312  // if(nameVolume == "C3EF" || nameVolume == "C4EF" || nameVolume == "C3HF" || nameVolume == "C4HF") {
313 
314  double bThreshold = 0.67;
315  double nMedium = 1.4925;
316  // double photEnSpectrDL = (1./400.nm-1./700.nm)*10000000.cm/nm; /* cm-1 */
317  // double photEnSpectrDL = 10714.285714;
318 
319  double photEnSpectrDE = 1.24; /* see below */
320  /* E = 2pi*(1./137.)*(eV*cm/370.)/lambda = */
321  /* = 12.389184*(eV*cm)/lambda */
322  /* Emax = 12.389184*(eV*cm)/400nm*10-7cm/nm = 3.01 eV */
323  /* Emin = 12.389184*(eV*cm)/700nm*10-7cm/nm = 1.77 eV */
324  /* delE = Emax - Emin = 1.24 eV --> */
325  /* */
326  /* default for Castor nameVolume == "CASF" or (C3TF & C4TF) */
327 
328  double thFullRefl = 23.; /* 23.dergee */
329  double thFullReflRad = thFullRefl*pi/180.;
330 
331  /* default for Castor nameVolume == "CASF" or (C3TF & C4TF) */
332  double thFibDir = 45.; /* .dergee */
333  /* for test HF geometry volumes:
334  if(nameVolume == "GF2Q" || nameVolume == "GFNQ" ||
335  nameVolume == "GR2Q" || nameVolume == "GRNQ")
336  thFibDir = 0.0; // .dergee
337  */
338  double thFibDirRad = thFibDir*pi/180.;
339  /* */
340  /* */
341 
342  // at which theta the point is located:
343  // double th1 = hitPoint.theta();
344 
345  // theta of charged particle in LabRF(hit momentum direction):
346  double costh =hit_mom.z()/sqrt(hit_mom.x()*hit_mom.x()+
347  hit_mom.y()*hit_mom.y()+
348  hit_mom.z()*hit_mom.z());
349  if (zint < 0) costh = -costh;
350  double th = acos(std::min(std::max(costh,double(-1.)),double(1.)));
351 
352  // just in case (can do bot use):
353  if (th < 0.) th += twopi;
354 
355 
356 
357  // theta of cone with Cherenkov photons w.r.t.direction of charged part.:
358  double costhcher =1./(nMedium*beta);
359  double thcher = acos(std::min(std::max(costhcher,double(-1.)),double(1.)));
360 
361  // diff thetas of charged part. and quartz direction in LabRF:
362  double DelFibPart = fabs(th - thFibDirRad);
363 
364  // define real distances:
365  double d = fabs(tan(th)-tan(thFibDirRad));
366 
367  // double a = fabs(tan(thFibDirRad)-tan(thFibDirRad+thFullReflRad));
368  // double r = fabs(tan(th)-tan(th+thcher));
369 
370  double a = tan(thFibDirRad)+tan(fabs(thFibDirRad-thFullReflRad));
371  double r = tan(th)+tan(fabs(th-thcher));
372 
373 
374  // define losses d_qz in cone of full reflection inside quartz direction
375  double d_qz;
376 #ifdef debugLog
377  double variant;
378 #endif
379  if(DelFibPart > (thFullReflRad + thcher) ) {
380  d_qz = 0.;
381 #ifdef debugLog
382  variant=0.;
383 #endif
384  }
385  // if(d > (r+a) ) {d_qz = 0.; variant=0.;}
386  else {
387  if((th + thcher) < (thFibDirRad+thFullReflRad) &&
388  (th - thcher) > (thFibDirRad-thFullReflRad)) {
389  d_qz = 1.;
390 #ifdef debugLog
391  variant=1.;
392 #endif
393  }
394  // if((DelFibPart + thcher) < thFullReflRad ) {d_qz = 1.; variant=1.;}
395  // if((d+r) < a ) {d_qz = 1.; variant=1.;}
396  else {
397  if((thFibDirRad + thFullReflRad) < (th + thcher) &&
398  (thFibDirRad - thFullReflRad) > (th - thcher) ) {
399  // if((thcher - DelFibPart ) > thFullReflRad )
400  // if((r-d ) > a )
401  d_qz = 0.;
402 #ifdef debugLog
403  variant=2.;
404 #endif
405  } else {
406  // if((thcher + DelFibPart ) > thFullReflRad &&
407  // thcher < (DelFibPart+thFullReflRad) )
408  // {
409 #ifdef debugLog
410  variant=3.;
411 #endif
412 
413  // use crossed length of circles(cone projection)
414  // dC1/dC2 :
415  double arg_arcos = 0.;
416  double tan_arcos = 2.*a*d;
417  if(tan_arcos != 0.) arg_arcos =(r*r-a*a-d*d)/tan_arcos;
418  arg_arcos = fabs(arg_arcos);
419  double th_arcos = acos(std::min(std::max(arg_arcos,double(-1.)),double(1.)));
420  d_qz = fabs(th_arcos/pi/2.);
421 
422  // }
423  // else
424  // {
425  // d_qz = 0.; variant=4.;
426  //#ifdef debugLog
427  // std::cout <<" ===============>variant 4 information: <===== " <<std::endl;
428  // std::cout <<" !!!!!!!!!!!!!!!!!!!!!! variant = " << variant <<std::endl;
429  //#endif
430  //
431  // }
432  }
433  }
434  }
435 
436 
437  // !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
438 
439  if(charge != 0. && beta > bThreshold ) {
440 
441  meanNCherPhot = 370.*charge*charge*
442  ( 1. - 1./(nMedium*nMedium*beta*beta) )*
443  photEnSpectrDE*stepl;
444 
445  const double scale = (isHad ? non_compensation_factor : 1.0);
446  G4int poissNCherPhot = (G4int) G4Poisson(meanNCherPhot * scale);
447 
448  if(poissNCherPhot < 0) poissNCherPhot = 0;
449 
450  double effPMTandTransport = 0.19;
451  double ReflPower = 0.1;
452  double proba = d_qz + (1-d_qz)*ReflPower;
453  NCherPhot = poissNCherPhot*effPMTandTransport*proba*0.307;
454 
455 
456 #ifdef debugLog
457  double thgrad = th*180./pi;
458  double thchergrad = thcher*180./pi;
459  double DelFibPartgrad = DelFibPart*180./pi;
460  LogDebug("ForwardSim") << " ==============================> start all "
461  << "information:<========= \n" << " =====> for "
462  << "test:<=== \n" << " variant = " << variant
463  << "\n thgrad = " << thgrad << "\n thchergrad "
464  << "= " << thchergrad << "\n DelFibPartgrad = "
465  << DelFibPartgrad << "\n d_qz = " << d_qz
466  << "\n =====> Start Step Information <=== \n"
467  << " ===> calo preStepPoint info <=== \n"
468  << " hitPoint = " << hitPoint << "\n"
469  << " hitMom = " << hit_mom << "\n"
470  << " stepControlFlag = " << stepControlFlag
471  // << "\n curprocess = " << curprocess << "\n"
472  // << " nameProcess = " << nameProcess
473  << "\n charge = " << charge << "\n"
474  << " beta = " << beta << "\n"
475  << " bThreshold = " << bThreshold << "\n"
476  << " thgrad =" << thgrad << "\n"
477  << " effPMTandTransport=" << effPMTandTransport
478  // << "\n volume = " << name
479  << "\n nameVolume = " << nameVolume << "\n"
480  << " nMedium = " << nMedium << "\n"
481  // << " rad length = " << rad << "\n"
482  // << " material = " << mat << "\n"
483  << " stepl = " << stepl << "\n"
484  << " photEnSpectrDE = " << photEnSpectrDE <<"\n"
485  << " edep = " << edep << "\n"
486  << " ===> calo theTrack info <=== " << "\n"
487  << " particleType = " << particleType << "\n"
488  << " primaryID = " << primaryID << "\n"
489  << " entot= " << theTrack->GetTotalEnergy() << "\n"
490  << " vert_eta= " << eta << "\n"
491  << " vert_phi= " << phi << "\n"
492  << " vert_mom= " << vert_mom << "\n"
493  << " ===> calo hit preStepPointinfo <=== "<<"\n"
494  << " local point = " << localPoint << "\n"
495  << " ==============================> final info"
496  << ": <=== \n"
497  << " meanNCherPhot = " << meanNCherPhot << "\n"
498  << " poissNCherPhot = " << poissNCherPhot <<"\n"
499  << " NCherPhot = " << NCherPhot;
500 #endif
501 
502  // Included by WC
503  // std::cout << "\n volume = " << name
504  // << "\n nameVolume = " << nameVolume << "\n"
505  // << "\n postvolume = " << postname
506  // << "\n postnameVolume = " << postnameVolume << "\n"
507  // << "\n particleType = " << particleType
508  // << "\n primaryID = " << primaryID << "\n";
509 
510  }
511  }
512 
513 
514 #ifdef debugLog
515  LogDebug("ForwardSim") << "CastorSD:: " << nameVolume
516  // << " Light Collection Efficiency " << weight
517  << " Weighted Energy Deposit " << edep/MeV
518  << " MeV\n";
519 #endif
520  // Temporary member for testing purpose only...
521  // unit_id = setDetUnitId(aStep);
522  // if(NCherPhot != 0) std::cout << "\n UnitID = " << unit_id << " ; NCherPhot = " << NCherPhot ;
523 
524  return NCherPhot;
525 }
526 
527 //=======================================================================================
528 
529 uint32_t CastorSD::setDetUnitId(const G4Step* aStep) {
530  return (numberingScheme == nullptr ? 0 : numberingScheme->getUnitID(aStep));
531 }
532 
533 //=======================================================================================
534 
536 
537  if (scheme != nullptr) {
538  edm::LogInfo("ForwardSim") << "CastorSD: updates numbering scheme for "
539  << GetName();
540  if (numberingScheme) delete numberingScheme;
541  numberingScheme = scheme;
542  }
543 }
544 
545 //=======================================================================================
546 
547 int CastorSD::setTrackID (G4Step* aStep) {
548 
549  theTrack = aStep->GetTrack();
550 
551  TrackInformation * trkInfo = (TrackInformation *)(theTrack->GetUserInformation());
552  int primaryID = trkInfo->getIDonCaloSurface();
553  if (primaryID == 0) {
554 #ifdef debugLog
555  edm::LogWarning("ForwardSim") << "CastorSD: Problem with primaryID **** set by force "
556  << "to TkID **** " << theTrack->GetTrackID();
557 #endif
558  primaryID = theTrack->GetTrackID();
559  }
560 
561  if (primaryID != previousID.trackID()) {
562  double etrack = preStepPoint->GetKineticEnergy();
563  resetForNewPrimary(preStepPoint->GetPosition(), etrack);
564  }
565 
566  return primaryID;
567 }
568 
569 //=======================================================================================
570 
571 uint32_t CastorSD::rotateUnitID(uint32_t unitID, G4Track* track, const CastorShowerEvent& shower) {
572 // ==============================================================
573 //
574 // o Exploit Castor phi symmetry to return newUnitID for
575 // shower hits based on track phi coordinate
576 //
577 // ==============================================================
578 
579  // Get 'track' phi:
580  double trackPhi = track->GetPosition().phi();
581  if(trackPhi<0) trackPhi += 2*M_PI ;
582  // Get phi from primary that gave rise to SL 'shower':
583  double showerPhi = shower.getPrimPhi();
584  if(showerPhi<0) showerPhi += 2*M_PI ;
585  // Delta phi:
586 
587  // Find the OctSector for which 'track' and 'shower' belong
588  int trackOctSector = (int) ( trackPhi / (M_PI/4) ) ;
589  int showerOctSector = (int) ( showerPhi / (M_PI/4) ) ;
590 
591  uint32_t newUnitID;
592  uint32_t sec = ( ( unitID>>4 ) & 0xF ) ;
593  uint32_t complement = ( unitID & 0xFFFFFF0F ) ;
594 
595  // edm::LogInfo("ForwardSim") << "\n CastorSD::rotateUnitID: "
596  // << "\n unitID = " << unitID
597  // << "\n sec = " << sec
598  // << "\n complement = " << complement ;
599 
600  // Get 'track' z:
601  double trackZ = track->GetPosition().z();
602 
603  int aux ;
604  int dSec = 2*(trackOctSector - showerOctSector) ;
605  // if(trackZ<0) // Good for revision 1.8 of CastorNumberingScheme
606  if(trackZ>0) // Good for revision 1.9 of CastorNumberingScheme
607  {
608  int sec1 = sec-dSec;
609  // sec -= dSec ;
610  if(sec1<0) sec1 += 16;
611  if(sec1>15) sec1 -= 16;
612  sec = (uint32_t)(sec1);
613  } else {
614  if( dSec<0 ) sec += 16 ;
615  sec += dSec ;
616  aux = (int) (sec/16) ;
617  sec -= aux*16 ;
618  }
619  sec = sec<<4 ;
620  newUnitID = complement | sec ;
621 
622 #ifdef debugLog
623  if(newUnitID != unitID) {
624  LogDebug("ForwardSim") << "\n CastorSD::rotateUnitID: "
625  << "\n unitID = " << unitID
626  << "\n newUnitID = " << newUnitID ;
627  }
628 #endif
629 
630  return newUnitID ;
631 
632 }
633 
634 //=======================================================================================
635 
636 void CastorSD::getFromLibrary (G4Step* aStep) {
637 
639 //
640 // Method to get hits from the Shower Library
641 //
642 // CastorShowerEvent hits returned by getShowerHits are used to
643 // replace the full simulation of the shower from theTrack
644 //
645 // "updateHit" save the Hits to a CaloG4Hit container
646 //
648 
649  preStepPoint = aStep->GetPreStepPoint();
650  theTrack = aStep->GetTrack();
651  bool ok;
652 
653  // **** Call method to retrieve hits from the ShowerLibrary ****
655 
656  double etrack = preStepPoint->GetKineticEnergy();
657  int primaryID = setTrackID(aStep);
658  // int primaryID = theTrack->GetTrackID();
659 
660  // Reset entry point for new primary
661  posGlobal = preStepPoint->GetPosition();
662  resetForNewPrimary(posGlobal, etrack);
663 
664  // Check whether track is EM or HAD
665  G4int particleCode = theTrack->GetDefinition()->GetPDGEncoding();
666  bool isEM , isHAD ;
667  if (particleCode==emPDG || particleCode==epPDG || particleCode==gammaPDG) {
668  isEM = true ; isHAD = false;
669  } else {
670  isEM = false; isHAD = true ;
671  }
672 
673 #ifdef debugLog
674  // edm::LogInfo("ForwardSim") << "\n CastorSD::getFromLibrary: "
675  LogDebug("ForwardSim") << "\n CastorSD::getFromLibrary: "
676  << hits.getNhit() << " hits for " << GetName() << " from "
677  << theTrack->GetDefinition()->GetParticleName() << " of "
678  << preStepPoint->GetKineticEnergy()/GeV << " GeV and trackID "
679  << theTrack->GetTrackID() ;
680 #endif
681 
682  // Scale to correct energy
683  double E_track = preStepPoint->GetTotalEnergy() ;
684  double E_SLhit = hits.getPrimE() * GeV ;
685  double scale = E_track/E_SLhit ;
686 
687  //Non compensation
688  if (isHAD){
689  scale=scale*non_compensation_factor; // if hadronic extend the scale with the non-compensation factor
690  } else {
691  scale=scale; // if electromagnetic, don't do anything
692  }
693 
694 
695 /* double theTrackEnergy = theTrack->GetTotalEnergy() ;
696 
697  if(fabs(theTrackEnergy-E_track)>10.) {
698  edm::LogInfo("ForwardSim") << "\n TrackID = " << theTrack->GetTrackID()
699  << "\n theTrackEnergy = " << theTrackEnergy
700  << "\n preStepPointEnergy = " << E_track ;
701  G4TrackVector tsec = *(aStep->GetSecondary());
702  for (unsigned int kk=0; kk<tsec.size(); kk++) {
703  edm::LogInfo("ForwardSim") << "CastorSD::getFromLibrary:"
704  << "\n tsec[" << kk << "]->GetTrackID() = "
705  << tsec[kk]->GetTrackID()
706  << " with energy "
707  << tsec[kk]->GetTotalEnergy() ;
708  }
709  }
710 */
711  // Loop over hits retrieved from the library
712  for (unsigned int i=0; i<hits.getNhit(); i++) {
713 
714  // Get nPhotoElectrons and set edepositEM / edepositHAD accordingly
715  double nPhotoElectrons = hits.getNphotons(i);
716  // Apply scaling
717  nPhotoElectrons *= scale ;
718  if(isEM) {
719  // edepositEM = nPhotoElectrons*GeV;
720  edepositEM = nPhotoElectrons;
721  edepositHAD = 0.;
722  } else if(isHAD) {
723  edepositEM = 0.;
724  edepositHAD = nPhotoElectrons;
725  // edepositHAD = nPhotoElectrons*GeV;
726  }
727 
728  // Get hit position and time
729  double time = hits.getTime(i);
730  // math::XYZPoint position = hits.getHitPosition(i);
731 
732  // Get hit detID
733  unsigned int unitID = hits.getDetID(i);
734 
735  // Make the detID "rotation" from one sector to another taking into account the
736  // sectors of the impinging particle (theTrack) and of the particle that produced
737  // the 'hits' retrieved from shower library
738  unsigned int rotatedUnitID = rotateUnitID(unitID , theTrack , hits);
739  currentID.setID(rotatedUnitID, time, primaryID, 0);
740  // currentID.setID(unitID, time, primaryID, 0);
741 
742  // check if it is in the same unit and timeslice as the previous one
743  if (currentID == previousID) {
745  } else {
746  if (!checkHit()) currentHit = createNewHit();
747  }
748  } // End of loop over hits
749 
750  //Now kill the current track
751  if (ok) {
752  theTrack->SetTrackStatus(fStopAndKill);
753 #ifdef debugLog
754  LogDebug("ForwardSim") << "CastorSD::getFromLibrary:"
755  << "\n \"theTrack\" with TrackID() = "
756  << theTrack->GetTrackID()
757  << " and with energy "
758  << theTrack->GetTotalEnergy()
759  << " has been set to be killed" ;
760 #endif
761  G4TrackVector tv = *(aStep->GetSecondary());
762  for (unsigned int kk=0; kk<tv.size(); kk++) {
763  if (tv[kk]->GetVolume() == preStepPoint->GetPhysicalVolume()) {
764  tv[kk]->SetTrackStatus(fStopAndKill);
765 #ifdef debugLog
766  LogDebug("ForwardSim") << "CastorSD::getFromLibrary:"
767  << "\n tv[" << kk << "]->GetTrackID() = "
768  << tv[kk]->GetTrackID()
769  << " with energy "
770  << tv[kk]->GetTotalEnergy()
771  << " has been set to be killed" ;
772 #endif
773  }
774  }
775  }
776 }
777 
#define LogDebug(id)
float edepositEM
Definition: CaloSD.h:122
const double beta
T getParameter(std::string const &) const
void initParticleTable(G4ParticleTable *)
G4int emPDG
Definition: CaloSD.h:137
int getCastorHitPID() const
const double GeV
Definition: MathUtil.h:16
void updateHit(CaloG4Hit *)
Definition: CaloSD.cc:414
int getIDonCaloSurface() const
Definition: CaloSD.h:42
double getEnergyDeposit(G4Step *) override
Definition: CastorSD.cc:94
uint32_t setDetUnitId(const G4Step *step) override
Definition: CastorSD.cc:529
Geom::Theta< T > theta() const
double non_compensation_factor
Definition: CastorSD.h:53
G4ThreeVector posGlobal
Definition: CaloSD.h:114
CastorNumberingScheme * numberingScheme
Definition: CastorSD.h:46
bool useShowerLibrary
Definition: CastorSD.h:51
uint32_t rotateUnitID(uint32_t, G4Track *, const CastorShowerEvent &)
Definition: CastorSD.cc:571
#define nullptr
type of data representation of DDCompactView
Definition: DDCompactView.h:90
void getFromLibrary(G4Step *)
Definition: CastorSD.cc:636
unsigned int getDetID(int i)
G4bool checkHit()
Definition: CaloSD.cc:305
static TrackerG4SimHitNumberingScheme & numberingScheme(const DDCompactView &cpv, const GeometricDet &det)
G4LogicalVolume * lvC4HF
Definition: CastorSD.h:48
const Double_t pi
void setCastorHitPID(const int pid)
const double MeV
void resetForNewPrimary(const G4ThreeVector &, double)
Definition: CaloSD.cc:428
float edepositHAD
Definition: CaloSD.h:122
T sqrt(T t)
Definition: SSEVec.h:18
int trackID() const
Definition: CaloHitID.h:25
G4int epPDG
Definition: CaloSD.h:137
float getTime(int i)
void initRun() override
Definition: CastorSD.cc:83
G4int gammaPDG
Definition: CaloSD.h:137
int setTrackID(G4Step *)
Definition: CastorSD.cc:547
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
float getNphotons(int i)
G4LogicalVolume * lvC3HF
Definition: CastorSD.h:48
T min(T a, T b)
Definition: MathUtil.h:58
CaloHitID previousID
Definition: CaloSD.h:118
CaloG4Hit * currentHit
Definition: CaloSD.h:129
bool hasCastorHit() const
G4Track * theTrack
Definition: CaloSD.h:119
double energyThresholdSL
Definition: CastorSD.h:52
#define M_PI
void setID(uint32_t unitID, double timeSlice, int trackID, uint16_t depth=0)
Definition: CaloHitID.cc:44
G4StepPoint * preStepPoint
Definition: CaloSD.h:121
G4LogicalVolume * lvC4EF
Definition: CastorSD.h:48
static const G4LogicalVolume * GetVolume(const std::string &name)
CaloHitID currentID
Definition: CaloSD.h:118
float getPrimPhi() const
CastorSD(const std::string &, const DDCompactView &, const SensitiveDetectorCatalog &clg, edm::ParameterSet const &, const SimTrackManager *)
Definition: CastorSD.cc:30
~CastorSD() override
Definition: CastorSD.cc:77
G4LogicalVolume * lvCAST
Definition: CastorSD.h:49
CastorShowerLibrary * showerLibrary
Definition: CastorSD.h:47
CastorShowerEvent getShowerHits(G4Step *, bool &)
void setNumberingScheme(CastorNumberingScheme *scheme)
Definition: CastorSD.cc:535
G4LogicalVolume * lvC3EF
Definition: CastorSD.h:48
T dot(const Basic3DVector &v) const
Scalar product, or "dot" product, with a vector of same type.
double a
Definition: hdecay.h:121
virtual uint32_t getUnitID(const G4Step *aStep) const
CaloG4Hit * createNewHit()
Definition: CaloSD.cc:337
unsigned int getNhit()