9 #include "G4LogicalVolume.hh" 10 #include "G4LogicalVolumeStore.hh" 11 #include "G4Poisson.hh" 14 #include "G4VProcess.hh" 23 #include "G4PhysicalConstants.hh" 24 #include "G4SystemOfUnits.hh" 32 :
CaloSD(name, cpv, clg, p, manager) {
44 edm::LogInfo(
"EcalSim") <<
"Constructing a DreamSD with name " << GetName() <<
"\n" 45 <<
"DreamSD:: Use of Birks law is set to " << useBirk
46 <<
" with three constants kB = " << birk1 <<
", C1 = " << birk2 <<
", C2 = " << birk3 <<
"\n" 47 <<
" Slope for Light yield is set to " << slopeLY <<
"\n" 48 <<
" Parameterization of Cherenkov is set to " << doCherenkov_
60 double edep = aStep->GetTotalEnergyDeposit() *
weight;
66 LogDebug(
"EcalSim") <<
"DreamSD:: " << aStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() <<
" Side " <<
side 67 <<
" Light Collection Efficiency " << weight <<
" Weighted Energy Deposit " << edep /
MeV 76 DimensionMap::const_iterator ite =
xtalLMap.begin();
77 const G4LogicalVolume *lv = (ite->first);
78 G4Material *material = lv->GetMaterial();
79 edm::LogInfo(
"EcalSim") <<
"DreamSD::initRun: Initializes for material " << material->GetName() <<
" in " 84 edm::LogWarning(
"EcalSim") <<
"Couldn't retrieve material properties table\n" 85 <<
" Material = " << material->GetName();
93 const G4VTouchable *touch = aStep->GetPreStepPoint()->GetTouchable();
94 uint32_t
id = (touch->GetReplicaNumber(1)) * 10 + (touch->GetReplicaNumber(0));
105 G4String attribute =
"ReadOutName";
110 const G4LogicalVolumeStore *lvs = G4LogicalVolumeStore::GetInstance();
111 std::vector<G4LogicalVolume *>::const_iterator lvcite;
117 G4LogicalVolume *lv =
nullptr;
118 for (lvcite = lvs->begin(); lvcite != lvs->end(); lvcite++)
119 if ((*lvcite)->GetName() ==
name) {
123 LogDebug(
"EcalSim") <<
"DreamSD::initMap (for " << sd <<
"): Solid " << name <<
" Shape " << sol.
shape()
124 <<
" Parameter 0 = " << paras[0] <<
" Logical Volume " << lv;
125 double length = 0,
width = 0;
128 length = 2.0 * paras.back();
129 width = 2.0 * paras.front();
133 LogDebug(
"EcalSim") <<
"DreamSD: Length Table for " << attribute <<
" = " << sd <<
":";
134 DimensionMap::const_iterator ite =
xtalLMap.begin();
136 for (; ite !=
xtalLMap.end(); ite++, i++) {
137 G4String
name =
"Unknown";
138 if (ite->first !=
nullptr)
139 name = (ite->first)->GetName();
140 LogDebug(
"EcalSim") <<
" " << i <<
" " << ite->first <<
" " << name <<
" L = " << ite->second.first
141 <<
" W = " << ite->second.second;
147 auto const stepPoint = aStep->GetPreStepPoint();
148 auto const lv = stepPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
149 G4String nameVolume = lv->GetName();
152 G4ThreeVector localPoint =
setToLocal(stepPoint->GetPosition(), stepPoint->GetTouchable());
154 double localz = localPoint.x();
155 double dapd = 0.5 * crlength - flag * localz;
156 if (dapd >= -0.1 || dapd <= crlength + 0.1) {
160 edm::LogWarning(
"EcalSim") <<
"DreamSD: light coll curve : wrong distance " 161 <<
"to APD " << dapd <<
" crlength = " << crlength <<
" crystal name = " << nameVolume
162 <<
" z of localPoint = " << localz <<
" take weight = " <<
weight;
164 LogDebug(
"EcalSim") <<
"DreamSD, light coll curve : " << dapd <<
" crlength = " << crlength
165 <<
" crystal name = " << nameVolume <<
" z of localPoint = " << localz
166 <<
" take weight = " <<
weight;
173 DimensionMap::const_iterator ite =
xtalLMap.find(lv);
175 length = ite->second.first;
182 DimensionMap::const_iterator ite =
xtalLMap.find(lv);
184 width = ite->second.second;
192 double cherenkovEnergy = 0;
194 return cherenkovEnergy;
195 G4Material *material = aStep->GetTrack()->GetMaterial();
199 if (Rindex ==
nullptr) {
201 return cherenkovEnergy;
206 int Rlength = Rindex->GetVectorLength() - 1;
207 double Pmin = Rindex->Energy(0);
208 double Pmax = Rindex->Energy(Rlength);
209 LogDebug(
"EcalSim") <<
"Material properties: " 211 <<
" Pmin = " << Pmin <<
" Pmax = " <<
Pmax;
214 const G4StepPoint *pPreStepPoint = aStep->GetPreStepPoint();
215 const G4StepPoint *pPostStepPoint = aStep->GetPostStepPoint();
216 const G4ThreeVector &x0 = pPreStepPoint->GetPosition();
217 G4ThreeVector p0 = aStep->GetDeltaPosition().unit();
218 const G4DynamicParticle *aParticle = aStep->GetTrack()->GetDynamicParticle();
219 const double charge = aParticle->GetDefinition()->GetPDGCharge();
221 double beta = 0.5 * (pPreStepPoint->GetBeta() + pPostStepPoint->GetBeta());
222 double BetaInverse = 1.0 /
beta;
224 LogDebug(
"EcalSim") <<
"Particle properties: " 226 <<
" charge = " << charge <<
" beta = " <<
beta;
230 if (meanNumberOfPhotons <= 0.0) {
231 LogDebug(
"EcalSim") <<
"Mean number of photons is zero: " << meanNumberOfPhotons <<
", stopping here";
232 return cherenkovEnergy;
236 meanNumberOfPhotons *= aStep->GetStepLength();
239 int numPhotons =
static_cast<int>(G4Poisson(meanNumberOfPhotons));
241 if (numPhotons <= 0) {
242 LogDebug(
"EcalSim") <<
"Poission number of photons is zero: " << numPhotons <<
", stopping here";
243 return cherenkovEnergy;
247 double dp = Pmax - Pmin;
248 double maxCos = BetaInverse / (*Rindex)[Rlength];
249 double maxSin2 = (1.0 - maxCos) * (1.0 + maxCos);
252 for (
int iPhoton = 0; iPhoton < numPhotons; ++iPhoton) {
255 double sampledMomentum, sampledRI;
256 double cosTheta, sin2Theta;
260 randomNumber = G4UniformRand();
261 sampledMomentum = Pmin + randomNumber * dp;
262 sampledRI = Rindex->Value(sampledMomentum);
263 cosTheta = BetaInverse / sampledRI;
265 sin2Theta = (1.0 - cosTheta) * (1.0 + cosTheta);
266 randomNumber = G4UniformRand();
268 }
while (randomNumber * maxSin2 > sin2Theta);
272 randomNumber = G4UniformRand();
274 double phi = twopi * randomNumber;
275 double sinPhi =
sin(phi);
276 double cosPhi =
cos(phi);
281 double sinTheta =
sqrt(sin2Theta);
282 double px = sinTheta * cosPhi;
283 double py = sinTheta * sinPhi;
284 double pz = cosTheta;
285 G4ThreeVector photonDirection(px, py, pz);
288 photonDirection.rotateUz(p0);
291 randomNumber = G4UniformRand();
292 G4ThreeVector photonPosition = x0 + randomNumber * aStep->GetDeltaPosition();
293 G4ThreeVector photonMomentum = sampledMomentum * photonDirection;
298 return cherenkovEnergy;
306 const G4Material *aMaterial,
307 const G4MaterialPropertyVector *Rindex) {
308 const double rFact = 369.81 / (eV * cm);
313 double BetaInverse = 1. /
beta;
320 int Rlength = Rindex->GetVectorLength() - 1;
321 double Pmin = Rindex->Energy(0);
322 double Pmax = Rindex->Energy(Rlength);
325 double nMin = (*Rindex)[0];
326 double nMax = (*Rindex)[Rlength];
331 double dp = 0., ge = 0., CAImin = 0.;
334 if (nMax < BetaInverse) {
338 else if (nMin > BetaInverse) {
348 Pmin = Rindex->Value(BetaInverse);
353 ge = CAImax - CAImin;
357 double numPhotons = rFact * charge / eplus * charge / eplus * (dp - ge * BetaInverse * BetaInverse);
359 LogDebug(
"EcalSim") <<
"@SUB=getAverageNumberOfPhotons" 360 <<
"CAImin = " << CAImin <<
"\n" 361 <<
"CAImax = " << CAImax <<
"\n" 362 <<
"dp = " << dp <<
", ge = " << ge <<
"\n" 363 <<
"numPhotons = " << numPhotons;
373 if (pbWO2Name != aMaterial->GetName()) {
375 <<
"expecting " << pbWO2Name <<
", got " << aMaterial->GetName();
379 G4MaterialPropertiesTable *table =
new G4MaterialPropertiesTable();
383 const int nEntries = 14;
384 double PhotonEnergy[nEntries] = {1.7712 * eV,
398 double RefractiveIndex[nEntries] = {2.17728,
413 table->AddProperty(
"RINDEX", PhotonEnergy, RefractiveIndex, nEntries);
414 aMaterial->SetMaterialPropertiesTable(table);
421 double currentRI = RefractiveIndex[
index];
422 double currentPM = PhotonEnergy[
index];
423 double currentCAI = 0.0;
425 double prevPM = currentPM;
426 double prevCAI = currentCAI;
427 double prevRI = currentRI;
428 while (++index < nEntries) {
429 currentRI = RefractiveIndex[
index];
430 currentPM = PhotonEnergy[
index];
431 currentCAI = 0.5 * (1.0 / (prevRI * prevRI) + 1.0 / (currentRI * currentRI));
432 currentCAI = prevCAI + (currentPM - prevPM) * currentCAI;
437 prevCAI = currentCAI;
441 LogDebug(
"EcalSim") <<
"Material properties set for " << aMaterial->GetName();
462 G4StepPoint *stepPoint = aStep->GetPreStepPoint();
463 G4LogicalVolume *lv = stepPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
464 G4String nameVolume = lv->GetName();
468 double dapd = 0.5 * crlength - x.x();
469 double y = p.y() / p.x() * dapd;
471 LogDebug(
"EcalSim") <<
"Distance to APD: " << dapd <<
" - y at APD: " <<
y;
478 double waveLength = p.mag() * 1.239e8;
482 LogDebug(
"EcalSim") <<
"Wavelength: " << waveLength <<
" - Energy: " << energy;
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
const std::vector< double > & parameters(void) const
Give the parameters of the solid.
const DDLogicalPart & logicalPart() const
The logical-part of the current node in the filtered-view.
double getPhotonEnergyDeposit_(const G4ParticleMomentum &p, const G4ThreeVector &x, const G4Step *aStep)
Returns energy deposit for a given photon.
Sin< T >::type sin(const T &t)
uint32_t setDetUnitId(const G4Step *) override
bool setPbWO2MaterialProperties_(G4Material *aMaterial)
Sets material properties at run-time...
const DDSolid & solid(void) const
Returns a reference object of the solid being the shape of this LogicalPart.
Compact representation of the geometrical detector hierarchy.
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
A DDSolid represents the shape of a part.
static double getEfficiency(const double &waveLengthNm)
Return efficiency for given photon wavelength (in nm)
bool next()
set current node to the next node in the filtered tree
double cherenkovDeposit_(const G4Step *aStep)
Returns the total energy due to Cherenkov radiation.
Cos< T >::type cos(const T &t)
Abs< T >::type abs(const T &t)
double getAverageNumberOfPhotons_(const double charge, const double beta, const G4Material *aMaterial, const G4MaterialPropertyVector *rIndex)
Returns average number of photons created by track.
DDSolidShape shape(void) const
The type of the solid.
double crystalWidth(G4LogicalVolume *) const
std::pair< double, double > Doubles
G4MaterialPropertiesTable * materialPropertiesTable
double curve_LY(const G4Step *, int)
G4ThreeVector setToLocal(const G4ThreeVector &, const G4VTouchable *) const
double crystalLength(G4LogicalVolume *) const
std::unique_ptr< G4PhysicsOrderedFreeVector > chAngleIntegrals_
Table of Cherenkov angle integrals vs photon momentum.
double getEnergyDeposit(const G4Step *) override
double getAttenuation(const G4Step *aStep, double birk1, double birk2, double birk3) const
DreamSD(const std::string &, const DDCompactView &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
bool firstChild()
set the current node to the first child ...
const std::string & name() const
Returns the name.
void initMap(const std::string &, const DDCompactView &)