CMS 3D CMS Logo

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