CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
DreamSD.cc
Go to the documentation of this file.
1 
8 
9 #include "G4LogicalVolumeStore.hh"
10 #include "G4LogicalVolume.hh"
11 #include "G4Step.hh"
12 #include "G4Track.hh"
13 #include "G4VProcess.hh"
14 #include "G4Poisson.hh"
15 
16 // Histogramming
19 #include <TTree.h>
20 
21 // Cherenkov
24 
25 //________________________________________________________________________________________
26 DreamSD::DreamSD(G4String name, const DDCompactView & cpv,
28  edm::ParameterSet const & p, const SimTrackManager* manager) :
29  CaloSD(name, cpv, clg, p, manager) {
30 
32  useBirk= m_EC.getParameter<bool>("UseBirkLaw");
33  doCherenkov_ = m_EC.getParameter<bool>("doCherenkov");
34  birk1 = m_EC.getParameter<double>("BirkC1")*(g/(MeV*cm2));
35  birk2 = m_EC.getParameter<double>("BirkC2");
36  birk3 = m_EC.getParameter<double>("BirkC3");
37  slopeLY= m_EC.getParameter<double>("SlopeLightYield");
38  readBothSide_ = m_EC.getUntrackedParameter<bool>("ReadBothSide", false);
39 
40  edm::LogInfo("EcalSim") << "Constructing a DreamSD with name " << GetName() << "\n"
41  << "DreamSD:: Use of Birks law is set to "
42  << useBirk << " with three constants kB = "
43  << birk1 << ", C1 = " << birk2 << ", C2 = "
44  << birk3 << "\n"
45  << " Slope for Light yield is set to "
46  << slopeLY << "\n"
47  << " Parameterization of Cherenkov is set to "
48  << doCherenkov_ << " and readout both sides is "
49  << readBothSide_;
50 
51  initMap(name,cpv);
52 
53  // Init histogramming
55 
56  if ( !tfile.isAvailable() )
57  throw cms::Exception("BadConfig") << "TFileService unavailable: "
58  << "please add it to config file";
59 
60  ntuple_ = tfile->make<TTree>("tree","Cherenkov photons");
61  if (doCherenkov_) {
62  ntuple_->Branch("nphotons",&nphotons_,"nphotons/I");
63  ntuple_->Branch("px",px_,"px[nphotons]/F");
64  ntuple_->Branch("py",py_,"py[nphotons]/F");
65  ntuple_->Branch("pz",pz_,"pz[nphotons]/F");
66  ntuple_->Branch("x",x_,"x[nphotons]/F");
67  ntuple_->Branch("y",y_,"y[nphotons]/F");
68  ntuple_->Branch("z",z_,"z[nphotons]/F");
69  }
70 
71 }
72 
73 //________________________________________________________________________________________
74 bool DreamSD::ProcessHits(G4Step * aStep, G4TouchableHistory *) {
75 
76  if (aStep == NULL) {
77  return true;
78  } else {
79  side = 1;
80  if (getStepInfo(aStep)) {
81  if (hitExists() == false && edepositEM+edepositHAD>0.)
83  if (readBothSide_) {
84  side = -1;
85  getStepInfo(aStep);
86  if (hitExists() == false && edepositEM+edepositHAD>0.)
88  }
89  }
90  }
91  return true;
92 }
93 
94 
95 //________________________________________________________________________________________
96 bool DreamSD::getStepInfo(G4Step* aStep) {
97 
98  preStepPoint = aStep->GetPreStepPoint();
99  theTrack = aStep->GetTrack();
100  G4String nameVolume = preStepPoint->GetPhysicalVolume()->GetName();
101 
102  // take into account light collection curve for crystals
103  double weight = 1.;
104  weight *= curve_LY(aStep, side);
105  if (useBirk) weight *= getAttenuation(aStep, birk1, birk2, birk3);
106  edepositEM = aStep->GetTotalEnergyDeposit() * weight;
107  LogDebug("EcalSim") << "DreamSD:: " << nameVolume << " Side " << side
108  <<" Light Collection Efficiency " << weight
109  << " Weighted Energy Deposit " << edepositEM/MeV
110  << " MeV";
111  // Get cherenkov contribution
112  if ( doCherenkov_ ) {
113  edepositHAD = cherenkovDeposit_( aStep );
114  } else {
115  edepositHAD = 0;
116  }
117 
118  double time = (aStep->GetPostStepPoint()->GetGlobalTime())/nanosecond;
119  unsigned int unitID= setDetUnitId(aStep);
120  if (side < 0) unitID++;
121  TrackInformation * trkInfo = (TrackInformation *)(theTrack->GetUserInformation());
122  int primaryID;
123 
124  if (trkInfo)
125  primaryID = trkInfo->getIDonCaloSurface();
126  else
127  primaryID = 0;
128 
129  if (primaryID == 0) {
130  edm::LogWarning("EcalSim") << "CaloSD: Problem with primaryID **** set by "
131  << "force to TkID **** "
132  << theTrack->GetTrackID() << " in Volume "
133  << preStepPoint->GetTouchable()->GetVolume(0)->GetName();
134  primaryID = theTrack->GetTrackID();
135  }
136 
137  bool flag = (unitID > 0);
138  G4TouchableHistory* touch =(G4TouchableHistory*)(theTrack->GetTouchable());
139  if (flag) {
140  currentID.setID(unitID, time, primaryID, 0);
141 
142  LogDebug("EcalSim") << "CaloSD:: GetStepInfo for"
143  << " PV " << touch->GetVolume(0)->GetName()
144  << " PVid = " << touch->GetReplicaNumber(0)
145  << " MVid = " << touch->GetReplicaNumber(1)
146  << " Unit " << currentID.unitID()
147  << " Edeposit = " << edepositEM << " " << edepositHAD;
148  } else {
149  LogDebug("EcalSim") << "CaloSD:: GetStepInfo for"
150  << " PV " << touch->GetVolume(0)->GetName()
151  << " PVid = " << touch->GetReplicaNumber(0)
152  << " MVid = " << touch->GetReplicaNumber(1)
153  << " Unit " << std::hex << unitID << std::dec
154  << " Edeposit = " << edepositEM << " " << edepositHAD;
155  }
156  return flag;
157 
158 }
159 
160 
161 //________________________________________________________________________________________
163 
164  // Get the material and set properties if needed
165  DimensionMap::const_iterator ite = xtalLMap.begin();
166  G4LogicalVolume* lv = (ite->first);
167  G4Material* material = lv->GetMaterial();
168  edm::LogInfo("EcalSim") << "DreamSD::initRun: Initializes for material "
169  << material->GetName() << " in " << lv->GetName();
170  materialPropertiesTable = material->GetMaterialPropertiesTable();
171  if ( !materialPropertiesTable ) {
172  if ( !setPbWO2MaterialProperties_( material ) ) {
173  edm::LogWarning("EcalSim") << "Couldn't retrieve material properties table\n"
174  << " Material = " << material->GetName();
175  }
176  materialPropertiesTable = material->GetMaterialPropertiesTable();
177  }
178 }
179 
180 
181 //________________________________________________________________________________________
182 uint32_t DreamSD::setDetUnitId(G4Step * aStep) {
183  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
184  uint32_t id = (touch->GetReplicaNumber(1))*10 + (touch->GetReplicaNumber(0));
185  LogDebug("EcalSim") << "DreamSD:: ID " << id;
186  return id;
187 }
188 
189 
190 //________________________________________________________________________________________
191 void DreamSD::initMap(G4String sd, const DDCompactView & cpv) {
192 
193  G4String attribute = "ReadOutName";
195  DDValue ddv(attribute,sd,0);
197  DDFilteredView fv(cpv);
198  fv.addFilter(filter);
199  fv.firstChild();
200 
201  const G4LogicalVolumeStore * lvs = G4LogicalVolumeStore::GetInstance();
202  std::vector<G4LogicalVolume *>::const_iterator lvcite;
203  bool dodet=true;
204  while (dodet) {
205  const DDSolid & sol = fv.logicalPart().solid();
206  std::vector<double> paras(sol.parameters());
207  G4String name = sol.name().name();
208  G4LogicalVolume* lv=0;
209  for (lvcite = lvs->begin(); lvcite != lvs->end(); lvcite++)
210  if ((*lvcite)->GetName() == name) {
211  lv = (*lvcite);
212  break;
213  }
214  LogDebug("EcalSim") << "DreamSD::initMap (for " << sd << "): Solid "
215  << name << " Shape " << sol.shape() <<" Parameter 0 = "
216  << paras[0] << " Logical Volume " << lv;
217  double length = 0, width = 0;
218  // Set length to be the largest size, width the smallest
219  std::sort( paras.begin(), paras.end() );
220  length = 2.0*paras.back();
221  width = 2.0*paras.front();
222  xtalLMap.insert( std::pair<G4LogicalVolume*,Doubles>(lv,Doubles(length,width)) );
223  dodet = fv.next();
224  }
225  LogDebug("EcalSim") << "DreamSD: Length Table for " << attribute << " = "
226  << sd << ":";
227  DimensionMap::const_iterator ite = xtalLMap.begin();
228  int i=0;
229  for (; ite != xtalLMap.end(); ite++, i++) {
230  G4String name = "Unknown";
231  if (ite->first != 0) name = (ite->first)->GetName();
232  LogDebug("EcalSim") << " " << i << " " << ite->first << " " << name
233  << " L = " << ite->second.first
234  << " W = " << ite->second.second;
235  }
236 }
237 
238 //________________________________________________________________________________________
239 double DreamSD::curve_LY(G4Step* aStep, int flag) {
240 
241  G4StepPoint* stepPoint = aStep->GetPreStepPoint();
242  G4LogicalVolume* lv = stepPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
243  G4String nameVolume= lv->GetName();
244 
245  double weight = 1.;
246  G4ThreeVector localPoint = setToLocal(stepPoint->GetPosition(),
247  stepPoint->GetTouchable());
248  double crlength = crystalLength(lv);
249  double localz = localPoint.x();
250  double dapd = 0.5 * crlength - flag*localz; // Distance from closest APD
251  if (dapd >= -0.1 || dapd <= crlength+0.1) {
252  if (dapd <= 100.)
253  weight = 1.0 + slopeLY - dapd * 0.01 * slopeLY;
254  } else {
255  edm::LogWarning("EcalSim") << "DreamSD: light coll curve : wrong distance "
256  << "to APD " << dapd << " crlength = "
257  << crlength << " crystal name = " << nameVolume
258  << " z of localPoint = " << localz
259  << " take weight = " << weight;
260  }
261  LogDebug("EcalSim") << "DreamSD, light coll curve : " << dapd
262  << " crlength = " << crlength
263  << " crystal name = " << nameVolume
264  << " z of localPoint = " << localz
265  << " take weight = " << weight;
266  return weight;
267 }
268 
269 //________________________________________________________________________________________
270 const double DreamSD::crystalLength(G4LogicalVolume* lv) const {
271 
272  double length= -1.;
273  DimensionMap::const_iterator ite = xtalLMap.find(lv);
274  if (ite != xtalLMap.end()) length = ite->second.first;
275  return length;
276 
277 }
278 
279 //________________________________________________________________________________________
280 const double DreamSD::crystalWidth(G4LogicalVolume* lv) const {
281 
282  double width= -1.;
283  DimensionMap::const_iterator ite = xtalLMap.find(lv);
284  if (ite != xtalLMap.end()) width = ite->second.second;
285  return width;
286 
287 }
288 
289 
290 //________________________________________________________________________________________
291 // Calculate total cherenkov deposit
292 // Inspired by Geant4's Cherenkov implementation
293 double DreamSD::cherenkovDeposit_( G4Step* aStep ) {
294 
295  double cherenkovEnergy = 0;
296  if (!materialPropertiesTable) return cherenkovEnergy;
297  G4Material* material = aStep->GetTrack()->GetMaterial();
298 
299  // Retrieve refractive index
300  const G4MaterialPropertyVector* Rindex = materialPropertiesTable->GetProperty("RINDEX");
301  if ( Rindex == NULL ) {
302  edm::LogWarning("EcalSim") << "Couldn't retrieve refractive index";
303  return cherenkovEnergy;
304  }
305 
306  LogDebug("EcalSim") << "Material properties: " << "\n"
307  << " Pmin = " << Rindex->GetMinPhotonEnergy()
308  << " Pmax = " << Rindex->GetMaxPhotonEnergy();
309 
310  // Get particle properties
311  G4StepPoint* pPreStepPoint = aStep->GetPreStepPoint();
312  G4StepPoint* pPostStepPoint = aStep->GetPostStepPoint();
313  G4ThreeVector x0 = pPreStepPoint->GetPosition();
314  G4ThreeVector p0 = aStep->GetDeltaPosition().unit();
315  const G4DynamicParticle* aParticle = aStep->GetTrack()->GetDynamicParticle();
316  const double charge = aParticle->GetDefinition()->GetPDGCharge();
317  // beta is averaged over step
318  double beta = 0.5*( pPreStepPoint->GetBeta() + pPostStepPoint->GetBeta() );
319  double BetaInverse = 1.0/beta;
320 
321  LogDebug("EcalSim") << "Particle properties: " << "\n"
322  << " charge = " << charge
323  << " beta = " << beta;
324 
325  // Now get number of photons generated in this step
326  double meanNumberOfPhotons = getAverageNumberOfPhotons_( charge, beta, material, Rindex );
327  if ( meanNumberOfPhotons <= 0.0 ) { // Don't do anything
328  LogDebug("EcalSim") << "Mean number of photons is zero: " << meanNumberOfPhotons
329  << ", stopping here";
330  return cherenkovEnergy;
331  }
332 
333  // number of photons is in unit of Geant4...
334  meanNumberOfPhotons *= aStep->GetStepLength();
335 
336  // Now get a poisson distribution
337  int numPhotons = static_cast<int>( G4Poisson(meanNumberOfPhotons) );
338  //edm::LogVerbatim("EcalSim") << "Number of photons = " << numPhotons;
339  if ( numPhotons <= 0 ) {
340  LogDebug("EcalSim") << "Poission number of photons is zero: " << numPhotons
341  << ", stopping here";
342  return cherenkovEnergy;
343  }
344 
345  // Material refraction properties
346  double Pmin = Rindex->GetMinPhotonEnergy();
347  double Pmax = Rindex->GetMaxPhotonEnergy();
348  double dp = Pmax - Pmin;
349  double maxCos = BetaInverse / Rindex->GetMaxProperty();
350  double maxSin2 = (1.0 - maxCos) * (1.0 + maxCos);
351 
352  // Finally: get contribution of each photon
353  for ( int iPhoton = 0; iPhoton<numPhotons; ++iPhoton ) {
354 
355  // Determine photon momentum
356  double randomNumber;
357  double sampledMomentum, sampledRI;
358  double cosTheta, sin2Theta;
359 
360  // sample a momentum (not sure why this is needed!)
361  do {
362  randomNumber = G4UniformRand();
363  sampledMomentum = Pmin + randomNumber * dp;
364  sampledRI = Rindex->GetProperty(sampledMomentum);
365  cosTheta = BetaInverse / sampledRI;
366 
367  sin2Theta = (1.0 - cosTheta)*(1.0 + cosTheta);
368  randomNumber = G4UniformRand();
369 
370  } while (randomNumber*maxSin2 > sin2Theta);
371 
372  // Generate random position of photon on cone surface
373  // defined by Theta
374  randomNumber = G4UniformRand();
375 
376  double phi = twopi*randomNumber;
377  double sinPhi = sin(phi);
378  double cosPhi = cos(phi);
379 
380  // Create photon momentum direction vector
381  // The momentum direction is still w.r.t. the coordinate system where the primary
382  // particle direction is aligned with the z axis
383  double sinTheta = sqrt(sin2Theta);
384  double px = sinTheta*cosPhi;
385  double py = sinTheta*sinPhi;
386  double pz = cosTheta;
387  G4ThreeVector photonDirection(px, py, pz);
388 
389  // Rotate momentum direction back to global (crystal) reference system
390  photonDirection.rotateUz(p0);
391 
392  // Create photon position and momentum
393  randomNumber = G4UniformRand();
394  G4ThreeVector photonPosition = x0 + randomNumber * aStep->GetDeltaPosition();
395  G4ThreeVector photonMomentum = sampledMomentum*photonDirection;
396 
397  // Collect energy on APD
398  cherenkovEnergy += getPhotonEnergyDeposit_( photonMomentum, photonPosition, aStep );
399 
400  // Ntuple variables
401  nphotons_ = numPhotons;
402  px_[iPhoton] = photonMomentum.x();
403  py_[iPhoton] = photonMomentum.y();
404  pz_[iPhoton] = photonMomentum.z();
405  x_[iPhoton] = photonPosition.x();
406  y_[iPhoton] = photonPosition.y();
407  z_[iPhoton] = photonPosition.z();
408  }
409 
410  // Fill ntuple
411  ntuple_->Fill();
412 
413 
414  return cherenkovEnergy;
415 
416 }
417 
418 
419 //________________________________________________________________________________________
420 // Returns number of photons produced per GEANT-unit (millimeter) in the current medium.
421 // From G4Cerenkov.cc
422 const double DreamSD::getAverageNumberOfPhotons_( const double charge,
423  const double beta,
424  const G4Material* aMaterial,
425  const G4MaterialPropertyVector* Rindex ) const
426 {
427  const G4double rFact = 369.81/(eV * cm);
428 
429  if( beta <= 0.0 ) return 0.0;
430 
431  double BetaInverse = 1./beta;
432 
433  // Vectors used in computation of Cerenkov Angle Integral:
434  // - Refraction Indices for the current material
435  // - new G4PhysicsOrderedFreeVector allocated to hold CAI's
436 
437  // Min and Max photon momenta
438  double Pmin = Rindex->GetMinPhotonEnergy();
439  double Pmax = Rindex->GetMaxPhotonEnergy();
440 
441  // Min and Max Refraction Indices
442  double nMin = Rindex->GetMinProperty();
443  double nMax = Rindex->GetMaxProperty();
444 
445  // Max Cerenkov Angle Integral
446  double CAImax = chAngleIntegrals_->GetMaxValue();
447 
448  double dp = 0., ge = 0., CAImin = 0.;
449 
450  // If n(Pmax) < 1/Beta -- no photons generated
451  if ( nMax < BetaInverse) { }
452 
453  // otherwise if n(Pmin) >= 1/Beta -- photons generated
454  else if (nMin > BetaInverse) {
455  dp = Pmax - Pmin;
456  ge = CAImax;
457  }
458  // If n(Pmin) < 1/Beta, and n(Pmax) >= 1/Beta, then
459  // we need to find a P such that the value of n(P) == 1/Beta.
460  // Interpolation is performed by the GetPhotonEnergy() and
461  // GetProperty() methods of the G4MaterialPropertiesTable and
462  // the GetValue() method of G4PhysicsVector.
463  else {
464  Pmin = Rindex->GetPhotonEnergy(BetaInverse);
465  dp = Pmax - Pmin;
466  // need boolean for current implementation of G4PhysicsVector
467  // ==> being phased out
468  bool isOutRange;
469  double CAImin = chAngleIntegrals_->GetValue(Pmin, isOutRange);
470  ge = CAImax - CAImin;
471 
472  }
473 
474  // Calculate number of photons
475  double numPhotons = rFact * charge/eplus * charge/eplus *
476  (dp - ge * BetaInverse*BetaInverse);
477 
478  LogDebug("EcalSim") << "@SUB=getAverageNumberOfPhotons"
479  << "CAImin = " << CAImin << "\n"
480  << "CAImax = " << CAImax << "\n"
481  << "dp = " << dp << ", ge = " << ge << "\n"
482  << "numPhotons = " << numPhotons;
483 
484 
485 
486  return numPhotons;
487 
488 }
489 
490 
491 //________________________________________________________________________________________
492 // Set lead tungstate material properties on the fly.
493 // Values from Ts42 detector construction
494 bool DreamSD::setPbWO2MaterialProperties_( G4Material* aMaterial ) {
495 
496  std::string pbWO2Name("E_PbWO4");
497  if ( pbWO2Name != aMaterial->GetName() ) { // Wrong material!
498  edm::LogWarning("EcalSim") << "This is not the right material: "
499  << "expecting " << pbWO2Name
500  << ", got " << aMaterial->GetName();
501  return false;
502  }
503 
504  G4MaterialPropertiesTable* table = new G4MaterialPropertiesTable();
505 
506  // Refractive index as a function of photon momentum
507  // FIXME: Should somehow put that in the configuration
508  const int nEntries = 14;
509  double PhotonEnergy[nEntries] = { 1.7712*eV, 1.8368*eV, 1.90745*eV, 1.98375*eV, 2.0664*eV,
510  2.15625*eV, 2.25426*eV, 2.3616*eV, 2.47968*eV, 2.61019*eV,
511  2.75521*eV, 2.91728*eV, 3.09961*eV, 3.30625*eV };
512  double RefractiveIndex[nEntries] = { 2.17728, 2.18025, 2.18357, 2.18753, 2.19285,
513  2.19813, 2.20441, 2.21337, 2.22328, 2.23619,
514  2.25203, 2.27381, 2.30282, 2.34666 };
515 
516  table->AddProperty( "RINDEX", PhotonEnergy, RefractiveIndex, nEntries );
517  aMaterial->SetMaterialPropertiesTable(table); // FIXME: could this leak? What does G4 do?
518 
519  // Calculate Cherenkov angle integrals:
520  // This is an ad-hoc solution (we hold it in the class, not in the material)
522  std::auto_ptr<G4PhysicsOrderedFreeVector>( new G4PhysicsOrderedFreeVector() );
523 
524  int index = 0;
525  double currentRI = RefractiveIndex[index];
526  double currentPM = PhotonEnergy[index];
527  double currentCAI = 0.0;
528  chAngleIntegrals_->InsertValues(currentPM, currentCAI);
529  double prevPM = currentPM;
530  double prevCAI = currentCAI;
531  double prevRI = currentRI;
532  while ( ++index < nEntries ) {
533  currentRI = RefractiveIndex[index];
534  currentPM = PhotonEnergy[index];
535  currentCAI = 0.5*(1.0/(prevRI*prevRI) + 1.0/(currentRI*currentRI));
536  currentCAI = prevCAI + (currentPM - prevPM) * currentCAI;
537 
538  chAngleIntegrals_->InsertValues(currentPM, currentCAI);
539 
540  prevPM = currentPM;
541  prevCAI = currentCAI;
542  prevRI = currentRI;
543  }
544 
545  LogDebug("EcalSim") << "Material properties set for " << aMaterial->GetName();
546 
547  return true;
548 
549 }
550 
551 
552 //________________________________________________________________________________________
553 // Calculate energy deposit of a photon on APD
554 // - simple tracing to APD position (straight line);
555 // - configurable reflection probability if not straight to APD;
556 // - APD response function
557 const double DreamSD::getPhotonEnergyDeposit_( const G4ThreeVector& p,
558  const G4ThreeVector& x,
559  const G4Step* aStep )
560 {
561 
562  double energy = 0;
563 
564  // Crystal dimensions
565 
566  //edm::LogVerbatim("EcalSim") << p << x;
567 
568  // 1. Check if this photon goes straight to the APD:
569  // - assume that APD is at x=xtalLength/2.0
570  // - extrapolate from x=x0 to x=xtalLength/2.0 using momentum in x-y
571 
572  G4StepPoint* stepPoint = aStep->GetPreStepPoint();
573  G4LogicalVolume* lv = stepPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
574  G4String nameVolume= lv->GetName();
575 
576  double crlength = crystalLength(lv);
577  double crwidth = crystalWidth(lv);
578  double dapd = 0.5 * crlength - x.x(); // Distance from closest APD
579  double y = p.y()/p.x()*dapd;
580 
581  LogDebug("EcalSim") << "Distance to APD: " << dapd
582  << " - y at APD: " << y;
583 
584  // Not straight: compute probability
585  if ( fabs(y)>crwidth*0.5 ) {
586 
587  }
588 
589  // 2. Retrieve efficiency for this wavelength (in nm, from MeV)
590  double waveLength = p.mag()*1.239e8;
591 
592 
593  energy = p.mag()*PMTResponse::getEfficiency(waveLength);
594 
595  LogDebug("EcalSim") << "Wavelength: " << waveLength << " - Energy: " << energy;
596 
597  return energy;
598 
599 }
G4ThreeVector setToLocal(G4ThreeVector, const G4VTouchable *)
Definition: CaloSD.cc:310
#define LogDebug(id)
float edepositEM
Definition: CaloSD.h:119
list table
Definition: asciidump.py:386
const double beta
const std::vector< double > & parameters() const
Don&#39;t use (only meant to be used by DDbox(), DDtub(), ...)
Definition: DDSolid.cc:153
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
long int flag
Definition: mlp_lapack.h:47
const DDLogicalPart & logicalPart() const
The logical-part of the current node in the filtered-view.
void addFilter(const DDFilter &, log_op op=AND)
const N & name() const
Definition: DDBase.h:88
bool readBothSide_
Definition: DreamSD.h:57
int getIDonCaloSurface() const
virtual G4bool getStepInfo(G4Step *aStep)
Definition: DreamSD.cc:96
std::pair< double, double > Doubles
Definition: DreamSD.h:36
Definition: CaloSD.h:42
DDSolidShape shape() const
The type of the solid.
Definition: DDSolid.cc:147
virtual bool ProcessHits(G4Step *step, G4TouchableHistory *tHistory)
Definition: DreamSD.cc:74
float y_[MAXPHOTONS]
Definition: DreamSD.h:71
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
#define NULL
Definition: scimark2.h:8
double birk1
Definition: DreamSD.h:58
bool setPbWO2MaterialProperties_(G4Material *aMaterial)
Sets material properties at run-time...
Definition: DreamSD.cc:494
std::auto_ptr< G4PhysicsOrderedFreeVector > chAngleIntegrals_
Table of Cherenkov angle integrals vs photon momentum.
Definition: DreamSD.h:65
double birk3
Definition: DreamSD.h:58
virtual uint32_t setDetUnitId(G4Step *)
Definition: DreamSD.cc:182
float px_[MAXPHOTONS]
Definition: DreamSD.h:70
type of data representation of DDCompactView
Definition: DDCompactView.h:81
TTree * ntuple_
Definition: DreamSD.h:68
double charge(const std::vector< uint8_t > &Ampls)
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e g
Definition: Activities.doc:4
A DDSolid represents the shape of a part.
Definition: DDSolid.h:42
virtual void initRun()
Definition: DreamSD.cc:162
const double getAverageNumberOfPhotons_(const double charge, const double beta, const G4Material *aMaterial, const G4MaterialPropertyVector *rIndex) const
Returns average number of photons created by track.
Definition: DreamSD.cc:422
float z_[MAXPHOTONS]
Definition: DreamSD.h:71
bool next()
set current node to the next node in the filtered tree
float edepositHAD
Definition: CaloSD.h:119
T sqrt(T t)
Definition: SSEVec.h:28
int nphotons_
Definition: DreamSD.h:69
float pz_[MAXPHOTONS]
Definition: DreamSD.h:70
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
int side
Definition: DreamSD.h:62
bool isAvailable() const
Definition: Service.h:47
DreamSD(G4String, const DDCompactView &, SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
Definition: DreamSD.cc:26
double birk2
Definition: DreamSD.h:58
double getAttenuation(G4Step *aStep, double birk1, double birk2, double birk3)
Definition: CaloSD.cc:480
CaloG4Hit * currentHit
Definition: CaloSD.h:126
G4MaterialPropertiesTable * materialPropertiesTable
Definition: DreamSD.h:66
G4Track * theTrack
Definition: CaloSD.h:116
const double getPhotonEnergyDeposit_(const G4ParticleMomentum &p, const G4ThreeVector &x, const G4Step *aStep)
Returns energy deposit for a given photon.
Definition: DreamSD.cc:557
const double crystalWidth(G4LogicalVolume *) const
Definition: DreamSD.cc:280
void setID(uint32_t unitID, double timeSlice, int trackID, uint16_t depth=0)
Definition: CaloHitID.cc:44
G4StepPoint * preStepPoint
Definition: CaloSD.h:118
double cherenkovDeposit_(G4Step *aStep)
Returns the total energy due to Cherenkov radiation.
Definition: DreamSD.cc:293
DimensionMap xtalLMap
Definition: DreamSD.h:60
tuple filter
USE THIS FOR SKIMMED TRACKS process.p = cms.Path(process.hltLevel1GTSeed*process.skimming*process.offlineBeamSpot*process.TrackRefitter2) OTHERWISE USE THIS.
Definition: align_tpl.py:86
CaloHitID currentID
Definition: CaloSD.h:115
float x_[MAXPHOTONS]
Definition: DreamSD.h:71
bool doCherenkov_
Definition: DreamSD.h:57
const DDSolid & solid() const
Returns a reference object of the solid being the shape of this LogicalPart.
double slopeLY
Definition: DreamSD.h:59
G4bool hitExists()
Definition: CaloSD.cc:319
bool firstChild()
set the current node to the first child ...
float py_[MAXPHOTONS]
Definition: DreamSD.h:70
T * make() const
make new ROOT object
void setCriteria(const DDValue &nameVal, comp_op, log_op l=AND, bool asString=true, bool merged=true)
Definition: DDFilter.cc:285
Definition: DDAxes.h:10
uint32_t unitID() const
Definition: CaloHitID.h:22
const double crystalLength(G4LogicalVolume *) const
Definition: DreamSD.cc:270
void initMap(G4String, const DDCompactView &)
Definition: DreamSD.cc:191
const std::string & name() const
Returns the name.
Definition: DDName.cc:87
double curve_LY(G4Step *, int)
Definition: DreamSD.cc:239
static const double getEfficiency(const double &waveLengthNm)
Return efficiency for given photon wavelength (in nm)
Definition: PMTResponse.cc:6
bool useBirk
Definition: DreamSD.h:57
CaloG4Hit * createNewHit()
Definition: CaloSD.cc:375
Definition: DDAxes.h:10
The DDGenericFilter is a runtime-parametrized Filter looking on DDSpecifcs.
Definition: DDFilter.h:37