11 #include "G4LogicalVolume.hh"
12 #include "G4LogicalVolumeStore.hh"
13 #include "G4Poisson.hh"
16 #include "G4VProcess.hh"
25 #include "G4PhysicalConstants.hh"
26 #include "G4SystemOfUnits.hh"
37 :
CaloSD(name, clg, p, manager), cpvDDD_(cpvDDD), cpvDD4hep_(cpvDD4hep) {
50 edm::LogVerbatim(
"EcalSim") <<
"Constructing a DreamSD with name " << GetName()
51 <<
"\nDreamSD:: Use of Birks law is set to " << useBirk_
52 <<
" with three constants kB = " << birk1_ <<
", C1 = " << birk2_ <<
", C2 = " << birk3_
53 <<
"\n Slope for Light yield is set to " << slopeLY_
54 <<
"\n Parameterization of Cherenkov is set to " << doCherenkov_
55 <<
", readout both sides is " << readBothSide_ <<
" and dd4hep flag " <<
dd4hep_;
58 const G4LogicalVolumeStore *lvs = G4LogicalVolumeStore::GetInstance();
60 for (
auto lvcite = lvs->begin(); lvcite != lvs->end(); ++lvcite, ++
k)
61 edm::LogVerbatim(
"EcalSim") <<
"Volume[" << k <<
"] " << (*lvcite)->GetName();
72 double edep = aStep->GetTotalEnergyDeposit() *
weight;
79 edm::LogVerbatim(
"EcalSim") <<
"DreamSD:: " << aStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() <<
" Side "
80 <<
side_ <<
" Light Collection Efficiency " << weight <<
" Weighted Energy Deposit "
89 DimensionMap::const_iterator ite =
xtalLMap_.begin();
90 const G4LogicalVolume *lv = (ite->first);
91 G4Material *material = lv->GetMaterial();
93 edm::LogVerbatim(
"EcalSim") <<
"DreamSD::initRun: Initializes for material " << material->GetName() <<
" in "
99 edm::LogWarning(
"EcalSim") <<
"Couldn't retrieve material properties table\n Material = " << material->GetName();
107 const G4VTouchable *touch = aStep->GetPreStepPoint()->GetTouchable();
108 uint32_t
id = (touch->GetReplicaNumber(1)) * 10 + (touch->GetReplicaNumber(0));
128 edm::LogVerbatim(
"EcalSim") <<
"DreamSD::initMap (for " << sd <<
"): Solid " << name <<
" Shape "
132 std::sort(paras.begin(), paras.end());
147 edm::LogVerbatim(
"EcalSim") <<
"DreamSD::initMap (for " << sd <<
"): Solid " << name <<
" Shape " << sol.
shape()
148 <<
" Parameter 0 = " << paras[0];
151 std::sort(paras.begin(), paras.end());
159 edm::LogVerbatim(
"EcalSim") <<
"DreamSD: Length Table for ReadOutName = " << sd <<
":";
161 DimensionMap::const_iterator ite =
xtalLMap_.begin();
163 for (; ite !=
xtalLMap_.end(); ite++, i++) {
164 G4String
name =
"Unknown";
165 if (ite->first !=
nullptr)
166 name = (ite->first)->GetName();
168 edm::LogVerbatim(
"EcalSim") <<
" " << i <<
" " << ite->first <<
" " << name <<
" L = " << ite->second.first
169 <<
" W = " << ite->second.second;
176 const G4LogicalVolumeStore *lvs = G4LogicalVolumeStore::GetInstance();
177 edm::LogVerbatim(
"EcalSim") <<
"LV Store with " << lvs->size() <<
" elements";
178 G4LogicalVolume *lv =
nullptr;
179 for (
auto lvcite = lvs->begin(); lvcite != lvs->end(); lvcite++) {
182 if (name == static_cast<std::string>(dd4hep::dd::noNamespace(namex))) {
188 edm::LogVerbatim(
"EcalSim") <<
"DreamSD::fillMap (for " << name <<
" Logical Volume " << lv <<
" Length " << length
189 <<
" Width " << width;
191 xtalLMap_.insert(std::pair<G4LogicalVolume *, Doubles>(lv,
Doubles(length, width)));
196 auto const stepPoint = aStep->GetPreStepPoint();
197 auto const lv = stepPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
198 G4String nameVolume = lv->GetName();
201 G4ThreeVector localPoint =
setToLocal(stepPoint->GetPosition(), stepPoint->GetTouchable());
203 double localz = localPoint.x();
204 double dapd = 0.5 * crlength - flag * localz;
205 if (dapd >= -0.1 || dapd <= crlength + 0.1) {
209 edm::LogWarning(
"EcalSim") <<
"DreamSD: light coll curve : wrong distance "
210 <<
"to APD " << dapd <<
" crlength = " << crlength <<
" crystal name = " << nameVolume
211 <<
" z of localPoint = " << localz <<
" take weight = " <<
weight;
214 edm::LogVerbatim(
"EcalSim") <<
"DreamSD, light coll curve : " << dapd <<
" crlength = " << crlength
215 <<
" crystal name = " << nameVolume <<
" z of localPoint = " << localz
216 <<
" take weight = " <<
weight;
224 DimensionMap::const_iterator ite =
xtalLMap_.find(lv);
226 length = ite->second.first;
233 DimensionMap::const_iterator ite =
xtalLMap_.find(lv);
235 width = ite->second.second;
243 double cherenkovEnergy = 0;
245 return cherenkovEnergy;
246 G4Material *material = aStep->GetTrack()->GetMaterial();
250 if (Rindex ==
nullptr) {
252 return cherenkovEnergy;
257 int Rlength = Rindex->GetVectorLength() - 1;
258 double Pmin = Rindex->Energy(0);
259 double Pmax = Rindex->Energy(Rlength);
261 edm::LogVerbatim(
"EcalSim") <<
"Material properties: \n Pmin = " << Pmin <<
" Pmax = " << Pmax;
264 const G4StepPoint *pPreStepPoint = aStep->GetPreStepPoint();
265 const G4StepPoint *pPostStepPoint = aStep->GetPostStepPoint();
266 const G4ThreeVector &x0 = pPreStepPoint->GetPosition();
267 G4ThreeVector p0 = aStep->GetDeltaPosition().unit();
268 const G4DynamicParticle *aParticle = aStep->GetTrack()->GetDynamicParticle();
269 const double charge = aParticle->GetDefinition()->GetPDGCharge();
271 double beta = 0.5 * (pPreStepPoint->GetBeta() + pPostStepPoint->GetBeta());
272 double BetaInverse = 1.0 /
beta;
275 edm::LogVerbatim(
"EcalSim") <<
"Particle properties: \n charge = " << charge <<
" beta = " <<
beta;
280 if (meanNumberOfPhotons <= 0.0) {
282 edm::LogVerbatim(
"EcalSim") <<
"Mean number of photons is zero: " << meanNumberOfPhotons <<
", stopping here";
284 return cherenkovEnergy;
288 meanNumberOfPhotons *= aStep->GetStepLength();
291 int numPhotons =
static_cast<int>(G4Poisson(meanNumberOfPhotons));
293 if (numPhotons <= 0) {
295 edm::LogVerbatim(
"EcalSim") <<
"Poission number of photons is zero: " << numPhotons <<
", stopping here";
297 return cherenkovEnergy;
301 double dp = Pmax - Pmin;
302 double maxCos = BetaInverse / (*Rindex)[Rlength];
303 double maxSin2 = (1.0 - maxCos) * (1.0 + maxCos);
306 for (
int iPhoton = 0; iPhoton < numPhotons; ++iPhoton) {
309 double sampledMomentum, sampledRI;
310 double cosTheta, sin2Theta;
314 randomNumber = G4UniformRand();
315 sampledMomentum = Pmin + randomNumber * dp;
316 sampledRI = Rindex->Value(sampledMomentum);
317 cosTheta = BetaInverse / sampledRI;
319 sin2Theta = (1.0 - cosTheta) * (1.0 + cosTheta);
320 randomNumber = G4UniformRand();
322 }
while (randomNumber * maxSin2 > sin2Theta);
326 randomNumber = G4UniformRand();
328 double phi = twopi * randomNumber;
335 double sinTheta =
sqrt(sin2Theta);
336 double px = sinTheta *
cosPhi;
337 double py = sinTheta *
sinPhi;
338 double pz = cosTheta;
339 G4ThreeVector photonDirection(px, py, pz);
342 photonDirection.rotateUz(p0);
345 randomNumber = G4UniformRand();
346 G4ThreeVector photonPosition = x0 + randomNumber * aStep->GetDeltaPosition();
347 G4ThreeVector photonMomentum = sampledMomentum * photonDirection;
352 return cherenkovEnergy;
360 const G4Material *aMaterial,
361 const G4MaterialPropertyVector *Rindex) {
362 const double rFact = 369.81 / (eV * cm);
367 double BetaInverse = 1. /
beta;
374 int Rlength = Rindex->GetVectorLength() - 1;
375 double Pmin = Rindex->Energy(0);
376 double Pmax = Rindex->Energy(Rlength);
379 double nMin = (*Rindex)[0];
380 double nMax = (*Rindex)[Rlength];
385 double dp = 0., ge = 0., CAImin = 0.;
388 if (nMax < BetaInverse) {
392 else if (nMin > BetaInverse) {
402 Pmin = Rindex->Value(BetaInverse);
407 ge = CAImax - CAImin;
411 double numPhotons = rFact * charge / eplus * charge / eplus * (dp - ge * BetaInverse * BetaInverse);
413 edm::LogVerbatim(
"EcalSim") <<
"@SUB=getAverageNumberOfPhotons\nCAImin = " << CAImin <<
"\nCAImax = " << CAImax
414 <<
"\ndp = " << dp <<
", ge = " << ge <<
"\nnumPhotons = " << numPhotons;
424 if (static_cast<std::string>(dd4hep::dd::noNamespace(name)) != pbWO2Name) {
426 <<
"expecting " << pbWO2Name <<
", got " << aMaterial->GetName();
430 G4MaterialPropertiesTable *
table =
new G4MaterialPropertiesTable();
434 const int nEntries = 14;
435 double PhotonEnergy[nEntries] = {1.7712 * eV,
449 double RefractiveIndex[nEntries] = {2.17728,
464 table->AddProperty(
"RINDEX", PhotonEnergy, RefractiveIndex, nEntries);
465 aMaterial->SetMaterialPropertiesTable(table);
472 double currentRI = RefractiveIndex[
index];
473 double currentPM = PhotonEnergy[
index];
474 double currentCAI = 0.0;
476 double prevPM = currentPM;
477 double prevCAI = currentCAI;
478 double prevRI = currentRI;
479 while (++index < nEntries) {
480 currentRI = RefractiveIndex[
index];
481 currentPM = PhotonEnergy[
index];
482 currentCAI = 0.5 * (1.0 / (prevRI * prevRI) + 1.0 / (currentRI * currentRI));
483 currentCAI = prevCAI + (currentPM - prevPM) * currentCAI;
488 prevCAI = currentCAI;
493 edm::LogVerbatim(
"EcalSim") <<
"Material properties set for " << aMaterial->GetName();
514 G4StepPoint *stepPoint = aStep->GetPreStepPoint();
515 G4LogicalVolume *lv = stepPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
516 G4String nameVolume = lv->GetName();
520 double dapd = 0.5 * crlength - x.x();
521 double y = p.y() / p.x() * dapd;
524 edm::LogVerbatim(
"EcalSim") <<
"Distance to APD: " << dapd <<
" - y at APD: " <<
y;
531 double waveLength = p.mag() * 1.239e8;
const DDCompactView * cpvDDD_
Log< level::Info, true > LogVerbatim
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.
static constexpr double k_ScaleFromDD4hepToG4
uint16_t *__restrict__ id
double getPhotonEnergyDeposit_(const G4ParticleMomentum &p, const G4ThreeVector &x, const G4Step *aStep)
Returns energy deposit for a given photon.
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
bool setPbWO2MaterialProperties_(G4Material *aMaterial)
Sets material properties at run-time...
static constexpr double k_ScaleFromDDDToG4
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.
const cms::DDSolidShape shape() 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.
uint16_t const *__restrict__ x
DDSolidShape shape(void) const
The type of the solid.
double crystalWidth(G4LogicalVolume *) const
const std::array< const cms::dd::NameValuePair< DDSolidShape >, 21 > DDSolidShapeMap
G4MaterialPropertiesTable * materialPropertiesTable_
std::string_view name() const
std::pair< double, double > Doubles
bool firstChild()
set the current node to the first child
double curve_LY(const G4Step *, int)
G4ThreeVector setToLocal(const G4ThreeVector &, const G4VTouchable *) const
DreamSD(const std::string &, const DDCompactView *, const cms::DDCompactView *, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
double crystalLength(G4LogicalVolume *) const
T getParameter(std::string const &) const
const cms::DDCompactView * cpvDD4hep_
double getEnergyDeposit(const G4Step *) override
double getAttenuation(const G4Step *aStep, double birk1, double birk2, double birk3) const
bool firstChild()
set the current node to the first child ...
void initMap(const std::string &)
const std::vector< double > parameters() const
extract shape parameters
Log< level::Warning, false > LogWarning
const std::string & name() const
Returns the name.