12 #include "G4LogicalVolume.hh" 13 #include "G4LogicalVolumeStore.hh" 14 #include "G4Poisson.hh" 17 #include "G4VProcess.hh" 26 #include "G4PhysicalConstants.hh" 27 #include "G4SystemOfUnits.hh" 38 :
CaloSD(
name, clg,
p, manager), cpvDDD_(cpvDDD), cpvDD4hep_(cpvDD4hep) {
47 dd4hep_ =
p.getParameter<
bool>(
"g4GeometryDD4hepSource");
51 edm::LogVerbatim(
"EcalSim") <<
"Constructing a DreamSD with name " << GetName()
52 <<
"\nDreamSD:: Use of Birks law is set to " <<
useBirk_ 53 <<
" with three constants kB = " <<
birk1_ <<
", C1 = " <<
birk2_ <<
", C2 = " <<
birk3_ 54 <<
"\n Slope for Light yield is set to " <<
slopeLY_ 55 <<
"\n Parameterization of Cherenkov is set to " <<
doCherenkov_ 59 const G4LogicalVolumeStore *lvs = G4LogicalVolumeStore::GetInstance();
61 for (
auto lvcite = lvs->begin(); lvcite != lvs->end(); ++lvcite, ++
k)
73 double edep = aStep->GetTotalEnergyDeposit() *
weight;
80 edm::LogVerbatim(
"EcalSim") <<
"DreamSD:: " << aStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() <<
" Side " 81 <<
side_ <<
" Light Collection Efficiency " <<
weight <<
" Weighted Energy Deposit " 82 << edep / CLHEP::MeV <<
" MeV";
90 DimensionMap::const_iterator ite =
xtalLMap_.begin();
91 const G4LogicalVolume *lv = (ite->first);
92 G4Material *material = lv->GetMaterial();
94 edm::LogVerbatim(
"EcalSim") <<
"DreamSD::initRun: Initializes for material " << material->GetName() <<
" in " 100 edm::LogWarning(
"EcalSim") <<
"Couldn't retrieve material properties table\n Material = " << material->GetName();
108 const G4VTouchable *touch = aStep->GetPreStepPoint()->GetTouchable();
109 uint32_t
id = (touch->GetReplicaNumber(1)) * 10 + (touch->GetReplicaNumber(0));
129 edm::LogVerbatim(
"EcalSim") <<
"DreamSD::initMap (for " << sd <<
"): Solid " <<
name <<
" Shape " 145 std::vector<double> paras(
sol.parameters());
148 edm::LogVerbatim(
"EcalSim") <<
"DreamSD::initMap (for " << sd <<
"): Solid " <<
name <<
" Shape " <<
sol.shape()
149 <<
" Parameter 0 = " << paras[0];
160 edm::LogVerbatim(
"EcalSim") <<
"DreamSD: Length Table for ReadOutName = " << sd <<
":";
162 DimensionMap::const_iterator ite =
xtalLMap_.begin();
165 G4String
name =
"Unknown";
166 if (ite->first !=
nullptr)
167 name = (ite->first)->GetName();
169 edm::LogVerbatim(
"EcalSim") <<
" " <<
i <<
" " << ite->first <<
" " <<
name <<
" L = " << ite->second.first
170 <<
" W = " << ite->second.second;
177 const G4LogicalVolumeStore *lvs = G4LogicalVolumeStore::GetInstance();
178 edm::LogVerbatim(
"EcalSim") <<
"LV Store with " << lvs->size() <<
" elements";
179 G4LogicalVolume *lv =
nullptr;
180 for (
auto lvcite = lvs->begin(); lvcite != lvs->end(); lvcite++) {
189 edm::LogVerbatim(
"EcalSim") <<
"DreamSD::fillMap (for " <<
name <<
" Logical Volume " << lv <<
" Length " << length
190 <<
" Width " <<
width;
197 auto const stepPoint = aStep->GetPreStepPoint();
198 auto const lv = stepPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
199 G4String nameVolume = lv->GetName();
202 G4ThreeVector localPoint =
setToLocal(stepPoint->GetPosition(), stepPoint->GetTouchable());
204 double localz = localPoint.x();
205 double dapd = 0.5 * crlength -
flag * localz;
206 if (dapd >= -0.1 || dapd <= crlength + 0.1) {
210 edm::LogWarning(
"EcalSim") <<
"DreamSD: light coll curve : wrong distance " 211 <<
"to APD " << dapd <<
" crlength = " << crlength <<
" crystal name = " << nameVolume
212 <<
" z of localPoint = " << localz <<
" take weight = " <<
weight;
215 edm::LogVerbatim(
"EcalSim") <<
"DreamSD, light coll curve : " << dapd <<
" crlength = " << crlength
216 <<
" crystal name = " << nameVolume <<
" z of localPoint = " << localz
217 <<
" take weight = " <<
weight;
225 DimensionMap::const_iterator ite =
xtalLMap_.find(lv);
227 length = ite->second.first;
234 DimensionMap::const_iterator ite =
xtalLMap_.find(lv);
236 width = ite->second.second;
244 double cherenkovEnergy = 0;
246 return cherenkovEnergy;
247 G4Material *material = aStep->GetTrack()->GetMaterial();
251 if (Rindex ==
nullptr) {
253 return cherenkovEnergy;
258 int Rlength = Rindex->GetVectorLength() - 1;
259 double Pmin = Rindex->Energy(0);
260 double Pmax = Rindex->Energy(Rlength);
265 const G4StepPoint *pPreStepPoint = aStep->GetPreStepPoint();
266 const G4StepPoint *pPostStepPoint = aStep->GetPostStepPoint();
267 const G4ThreeVector &x0 = pPreStepPoint->GetPosition();
268 G4ThreeVector p0 = aStep->GetDeltaPosition().unit();
269 const G4DynamicParticle *aParticle = aStep->GetTrack()->GetDynamicParticle();
270 const double charge = aParticle->GetDefinition()->GetPDGCharge();
272 double beta = 0.5 * (pPreStepPoint->GetBeta() + pPostStepPoint->GetBeta());
273 double BetaInverse = 1.0 /
beta;
281 if (meanNumberOfPhotons <= 0.0) {
283 edm::LogVerbatim(
"EcalSim") <<
"Mean number of photons is zero: " << meanNumberOfPhotons <<
", stopping here";
285 return cherenkovEnergy;
289 meanNumberOfPhotons *= aStep->GetStepLength();
292 int numPhotons =
static_cast<int>(G4Poisson(meanNumberOfPhotons));
294 if (numPhotons <= 0) {
296 edm::LogVerbatim(
"EcalSim") <<
"Poission number of photons is zero: " << numPhotons <<
", stopping here";
298 return cherenkovEnergy;
303 double maxCos = BetaInverse / (*Rindex)[Rlength];
304 double maxSin2 = (1.0 - maxCos) * (1.0 + maxCos);
307 for (
int iPhoton = 0; iPhoton < numPhotons; ++iPhoton) {
310 double sampledMomentum, sampledRI;
311 double cosTheta, sin2Theta;
315 randomNumber = G4UniformRand();
316 sampledMomentum = Pmin + randomNumber *
dp;
317 sampledRI = Rindex->Value(sampledMomentum);
318 cosTheta = BetaInverse / sampledRI;
320 sin2Theta = (1.0 - cosTheta) * (1.0 + cosTheta);
321 randomNumber = G4UniformRand();
323 }
while (randomNumber * maxSin2 > sin2Theta);
327 randomNumber = G4UniformRand();
329 double phi = twopi * randomNumber;
336 double sinTheta =
sqrt(sin2Theta);
339 double pz = cosTheta;
340 G4ThreeVector photonDirection(
px,
py, pz);
343 photonDirection.rotateUz(p0);
346 randomNumber = G4UniformRand();
347 G4ThreeVector photonPosition = x0 + randomNumber * aStep->GetDeltaPosition();
348 G4ThreeVector photonMomentum = sampledMomentum * photonDirection;
353 return cherenkovEnergy;
361 const G4Material *aMaterial,
362 const G4MaterialPropertyVector *Rindex) {
363 const double rFact = 369.81 / (eV * cm);
368 double BetaInverse = 1. /
beta;
375 int Rlength = Rindex->GetVectorLength() - 1;
376 double Pmin = Rindex->Energy(0);
377 double Pmax = Rindex->Energy(Rlength);
380 double nMin = (*Rindex)[0];
381 double nMax = (*Rindex)[Rlength];
386 double dp = 0., ge = 0., CAImin = 0.;
389 if (
nMax < BetaInverse) {
393 else if (
nMin > BetaInverse) {
403 Pmin = Rindex->Value(BetaInverse);
408 ge = CAImax - CAImin;
412 double numPhotons = rFact *
charge / eplus *
charge / eplus * (
dp - ge * BetaInverse * BetaInverse);
414 edm::LogVerbatim(
"EcalSim") <<
"@SUB=getAverageNumberOfPhotons\nCAImin = " << CAImin <<
"\nCAImax = " << CAImax
415 <<
"\ndp = " <<
dp <<
", ge = " << ge <<
"\nnumPhotons = " << numPhotons;
427 <<
"expecting " << pbWO2Name <<
", got " << aMaterial->GetName();
431 G4MaterialPropertiesTable *
table =
new G4MaterialPropertiesTable();
435 const int nEntries = 14;
436 double PhotonEnergy[nEntries] = {1.7712 * eV,
450 double RefractiveIndex[nEntries] = {2.17728,
465 table->AddProperty(
"RINDEX", PhotonEnergy, RefractiveIndex, nEntries);
466 aMaterial->SetMaterialPropertiesTable(
table);
473 double currentRI = RefractiveIndex[
index];
474 double currentPM = PhotonEnergy[
index];
475 double currentCAI = 0.0;
477 double prevPM = currentPM;
478 double prevCAI = currentCAI;
479 double prevRI = currentRI;
480 while (++
index < nEntries) {
481 currentRI = RefractiveIndex[
index];
482 currentPM = PhotonEnergy[
index];
483 currentCAI = 0.5 * (1.0 / (prevRI * prevRI) + 1.0 / (currentRI * currentRI));
484 currentCAI = prevCAI + (currentPM - prevPM) * currentCAI;
489 prevCAI = currentCAI;
494 edm::LogVerbatim(
"EcalSim") <<
"Material properties set for " << aMaterial->GetName();
515 G4StepPoint *stepPoint = aStep->GetPreStepPoint();
516 G4LogicalVolume *lv = stepPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
517 G4String nameVolume = lv->GetName();
521 double dapd = 0.5 * crlength -
x.x();
522 double y =
p.y() /
p.x() * dapd;
525 edm::LogVerbatim(
"EcalSim") <<
"Distance to APD: " << dapd <<
" - y at APD: " <<
y;
532 double waveLength =
p.mag() * 1.239e8;
const DDCompactView * cpvDDD_
Log< level::Info, true > LogVerbatim
double crystalLength(G4LogicalVolume *) const
T getParameter(std::string const &) const
static constexpr double k_ScaleFromDD4hepToG4
double getPhotonEnergyDeposit_(const G4ParticleMomentum &p, const G4ThreeVector &x, const G4Step *aStep)
Returns energy deposit for a given photon.
const cms::DDSolidShape shape() const
std::string name(Mapping a, V value)
Sin< T >::type sin(const T &t)
void fillMap(const std::string &, double, double)
uint32_t setDetUnitId(const G4Step *) override
G4ThreeVector setToLocal(const G4ThreeVector &, const G4VTouchable *) const
bool setPbWO2MaterialProperties_(G4Material *aMaterial)
Sets material properties at run-time...
static constexpr double k_ScaleFromDDDToG4
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
double crystalWidth(G4LogicalVolume *) const
A DDSolid represents the shape of a part.
T getUntrackedParameter(std::string const &, T const &) const
std::string_view name() const
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)
std::unique_ptr< G4PhysicsFreeVector > chAngleIntegrals_
Table of Cherenkov angle integrals vs photon momentum.
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.
const std::array< const cms::dd::NameValuePair< DDSolidShape >, 21 > DDSolidShapeMap
G4MaterialPropertiesTable * materialPropertiesTable_
std::pair< double, double > Doubles
bool firstChild()
set the current node to the first child
double curve_LY(const G4Step *, int)
const DDLogicalPart & logicalPart() const
The logical-part of the current node in the filtered-view.
DreamSD(const std::string &, const DDCompactView *, const cms::DDCompactView *, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
double getAttenuation(const G4Step *aStep, double birk1, double birk2, double birk3) const
const cms::DDCompactView * cpvDD4hep_
const DDSolid & solid(void) const
Returns a reference object of the solid being the shape of this LogicalPart.
double getEnergyDeposit(const G4Step *) override
bool firstChild()
set the current node to the first child ...
void initMap(const std::string &)
std::string noNameSpace(const std::string &name)
Log< level::Warning, false > LogWarning
const std::vector< double > parameters() const
extract shape parameters