21 #include "G4LogicalVolumeStore.hh"
22 #include "G4LogicalVolume.hh"
25 #include "G4VProcess.hh"
27 #include "G4SystemOfUnits.hh"
34 bool any(
const std::vector<T> &
v,
const T &what) {
35 return std::find(v.begin(), v.end(), what) != v.end();
41 CaloSD(name, cpv, clg, p, manager,
42 (float)(p.getParameter<edm::
ParameterSet>(
"ECalSD").getParameter<double>(
"TimeSliceUnit")),
43 p.getParameter<edm::
ParameterSet>(
"ECalSD").getParameter<bool>(
"IgnoreTrackID")),
88 std::vector<double> tempD =
getDDDArray(
"EnergyWeight",sv);
89 if (tempD.size() > 0) {
if (tempD[0] < 0.1) useWeight =
false; }
98 if (nullNS) scheme = 0;
99 else if (name ==
"EcalHitsEB")
104 else if (name ==
"EcalHitsEE")
109 else if (name ==
"EcalHitsES")
115 else {
edm::LogWarning(
"EcalSim") <<
"ECalSD: ReadoutName not supported\n";}
120 LogDebug(
"EcalSim") <<
"Constructing a ECalSD with name " << GetName();
123 edm::LogInfo(
"EcalSim") <<
"ECalSD:: Use of Birks law is set to "
124 << useBirk <<
" with three constants kB = "
125 << birk1 <<
", C1 = " << birk2 <<
", C2 = "
126 << birk3 <<
"\n Use of L3 parametrization "
127 << useBirkL3 <<
" with slope " << birkSlope
128 <<
" and cut off " << birkCut <<
"\n"
129 <<
" Slope for Light yield is set to "
132 edm::LogInfo(
"EcalSim") <<
"ECalSD:: energy deposit is not corrected "
133 <<
" by Birk or light yield curve";
137 <<
"\tprotons below " <<
kmaxProton <<
" MeV,"
139 <<
"\tions below " <<
kmaxIon <<
" MeV"
142 <<
"\n\tstoreRL" << storeRL
143 <<
"\tstoreLayerTimeSim " << storeLayerTimeSim
145 if (useWeight)
initMap(name,cpv);
158 G4Track*
theTrack = aStep->GetTrack();
159 double wt2 = theTrack->GetWeight();
160 G4String nameVolume =
preStepPoint->GetPhysicalVolume()->GetName();
167 int pdg = theTrack->GetDefinition()->GetPDGEncoding();
169 double ke = theTrack->GetKineticEnergy()/
MeV;
170 if (((pdg/1000000000 == 1 && ((pdg/10000)%100) > 0 &&
171 ((pdg/10)%100) > 0)) && (ke<
kmaxIon)) weight = 0;
172 if ((pdg == 2212) && (ke <
kmaxProton)) weight = 0;
173 if ((pdg == 2112) && (ke <
kmaxNeutron)) weight = 0;
176 LogDebug(
"EcalSim") <<
"Ignore Track " << theTrack->GetTrackID()
177 <<
" Type " << theTrack->GetDefinition()->GetParticleName()
178 <<
" Kinetic Energy " << ke <<
" MeV";
183 G4LogicalVolume* lv = aStep->GetPreStepPoint()->GetTouchable()->GetVolume(0)->GetLogicalVolume();
192 double edep = aStep->GetTotalEnergyDeposit()*weight*wt1;
209 if(wt2 > 0.0) { edep *= wt2; }
211 LogDebug(
"EcalSim") <<
"ECalSD:: " << nameVolume
212 <<
" Light Collection Efficiency " <<weight <<
":" <<wt1
213 <<
" Weighted Energy Deposit " << edep/
MeV <<
" MeV";
224 G4LogicalVolume* lv =
preStepPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
233 primaryID = aTrack->GetTrackID();
241 G4LogicalVolume* lv = aStep->GetPreStepPoint()->GetTouchable()->GetVolume(0)->GetLogicalVolume();
253 G4StepPoint* hitPoint = aStep->GetPreStepPoint();
254 G4LogicalVolume* lv = hitPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
257 G4ThreeVector localPoint =
setToLocal(hitPoint->GetPosition(),
258 hitPoint->GetTouchable());
260 double radl = hitPoint->GetMaterial()->GetRadlen();
261 double detz = (float)(0.5*crlength + localPoint.z());
262 thisX0 = (uint16_t)floor(detz/radl);
271 float layerSize = 1*cm;
275 if (aStep !=
NULL ) {
276 G4StepPoint* hitPoint = aStep->GetPreStepPoint();
277 G4LogicalVolume* lv = hitPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
278 G4ThreeVector localPoint =
setToLocal(hitPoint->GetPosition(),
279 hitPoint->GetTouchable());
283 const auto&
name = lv->GetName();
288 detz = (float)(0.5*crlength + localPoint.z());
290 detz = (float)(0.5*crlength - localPoint.z());
295 detz = (float)(0.5*crlength - localPoint.z());
297 detz = (float)(0.5*crlength + localPoint.z());
301 return 100+(int)detz/layerSize;
317 edm::LogInfo(
"EcalSim") <<
"EcalSD: updates numbering scheme for "
318 << GetName() <<
"\n";
327 G4String attribute =
"ReadOutName";
335 std::vector<G4LogicalVolume*> lvused;
336 const G4LogicalVolumeStore * lvs = G4LogicalVolumeStore::GetInstance();
337 std::vector<G4LogicalVolume *>::const_iterator lvcite;
338 std::map<std::string, G4LogicalVolume *> nameMap;
339 for (
auto lvi = lvs->begin(), lve = lvs->end(); lvi != lve; ++lvi)
340 nameMap.insert(std::make_pair((*lvi)->GetName(), *lvi));
346 G4LogicalVolume* lv = nameMap[lvname];
348 if (strncmp(lvname.c_str(),
depth1Name.c_str(), 4) == 0) {
352 LogDebug(
"EcalSim") <<
"ECalSD::initMap Logical Volume " << lvname
353 <<
" in Depth 1 volume list";
356 G4LogicalVolume* lvr = nameMap[lvname +
"_refl"];
360 LogDebug(
"EcalSim") <<
"ECalSD::initMap Logical Volume " << lvname <<
"_refl"
361 <<
" in Depth 1 volume list";
367 if (strncmp(lvname.c_str(),
depth2Name.c_str(), 4) == 0) {
371 LogDebug(
"EcalSim") <<
"ECalSD::initMap Logical Volume " << lvname
372 <<
" in Depth 2 volume list";
375 G4LogicalVolume* lvr = nameMap[lvname +
"_refl"];
379 LogDebug(
"EcalSim") <<
"ECalSD::initMap Logical Volume " << lvname <<
"_refl"
380 <<
" in Depth 2 volume list";
387 if (!
any(lvused,lv)) {
388 lvused.push_back(lv);
390 const std::vector<double> & paras = sol.
parameters();
392 LogDebug(
"EcalSim") <<
"ECalSD::initMap (for " << sd <<
"): Solid "
393 << lvname <<
" Shape " << sol.
shape()
394 <<
" Parameter 0 = "<< paras[0]
395 <<
" Logical Volume " << lv;
398 double dz = 2*paras[0];
399 xtalLMap.insert(std::pair<G4LogicalVolume*,double>(lv,dz));
400 lv = nameMap[lvname +
"_refl"];
402 xtalLMap.insert(std::pair<G4LogicalVolume*,double>(lv,dz));
409 LogDebug(
"EcalSim") <<
"ECalSD::initMap Logical Volume " << lvname
410 <<
" Material " << matname <<
" in noWeight list";
413 lv = nameMap[lvname];
417 LogDebug(
"EcalSim") <<
"ECalSD::initMap Logical Volume " << lvname
418 <<
" Material " << matname <<
" in noWeight list";
426 LogDebug(
"EcalSim") <<
"ECalSD: Length Table for " << attribute <<
" = "
428 std::map<G4LogicalVolume*,double>::const_iterator ite =
xtalLMap.begin();
430 for (; ite !=
xtalLMap.end(); ite++, i++) {
431 G4String
name =
"Unknown";
432 if (ite->first != 0) name = (ite->first)->GetName();
433 LogDebug(
"EcalSim") <<
" " << i <<
" " << ite->first <<
" " << name
434 <<
" L = " << ite->second;
441 G4StepPoint* stepPoint = aStep->GetPreStepPoint();
442 G4LogicalVolume* lv = stepPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
445 G4ThreeVector localPoint =
setToLocal(stepPoint->GetPosition(),
446 stepPoint->GetTouchable());
452 double depth = 0.5 * crlength + localPoint.z();
454 if (depth >= -0.1 || depth <= crlength+0.1)
458 double dapd = 0.5 * crlength - localPoint.z();
459 if (dapd >= -0.1 || dapd <= crlength+0.1) {
463 edm::LogWarning(
"EcalSim") <<
"ECalSD: light coll curve : wrong distance "
464 <<
"to APD " << dapd <<
" crlength = "
465 << crlength <<
" crystal name = " <<lv->GetName()
466 <<
" z of localPoint = " << localPoint.z()
467 <<
" take weight = " <<
weight;
476 std::map<G4LogicalVolume*,double>::const_iterator ite =
xtalLMap.find(lv);
477 if (ite !=
xtalLMap.end()) length = ite->second;
484 const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
485 int theSize = touch->GetHistoryDepth()+1;
489 for (
int ii = 0;
ii < theSize ;
ii++) {
492 LogDebug(
"EcalSim") <<
"ECalSD::getBaseNumber(): Adding level " <<
ii
493 <<
": " << touch->GetVolume(
ii)->GetName() <<
"["
494 << touch->GetReplicaNumber(
ii) <<
"]";
504 double charge = aStep->GetPreStepPoint()->GetCharge();
506 if (charge != 0. && aStep->GetStepLength() > 0) {
507 G4Material* mat = aStep->GetPreStepPoint()->GetMaterial();
508 double density = mat->GetDensity();
509 double dedx = aStep->GetTotalEnergyDeposit()/aStep->GetStepLength();
510 double rkb =
birk1/density;
514 else if (weight > 1.) weight = 1.;
517 LogDebug(
"EcalSim") <<
"ECalSD::getBirkL3 in " << mat->GetName()
518 <<
" Charge " << charge <<
" dE/dx " << dedx
519 <<
" Birk Const " << rkb <<
" Weight = " << weight
520 <<
" dE " << aStep->GetTotalEnergyDeposit();
531 LogDebug(
"EcalSim") <<
"ECalSD:getDDDArray called for " << str;
538 const std::vector<double> & fvec = value.
doubles();
541 std::vector<double> fvec;
550 LogDebug(
"EcalSim") <<
"ECalSD:getStringArray called for " << str;
557 const std::vector<std::string> & fvec = value.
strings();
560 std::vector<std::string> fvec;
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.
const std::vector< double > & doubles() const
a reference to the double-valued values stored in the given instance of DDValue
std::vector< std::string > getStringArray(const std::string &, const DDsvalues_type &)
std::vector< G4LogicalVolume * > useDepth1
void setLumies(double x, double y)
void getBaseNumber(const G4Step *)
void addFilter(const DDFilter &, DDLogOp op=DDLogOp::AND)
virtual double getEnergyDeposit(G4Step *)
bool any(const std::vector< T > &v, const T &what)
virtual int getTrackID(G4Track *)
const DDSolid & solid(void) const
Returns a reference object of the solid being the shape of this LogicalPart.
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
type of data representation of DDCompactView
bool DDfetch(const DDsvalues_type *, DDValue &)
helper for retrieving DDValues from DDsvalues_type *.
std::vector< G4LogicalVolume * > noWeight
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
std::vector< G4LogicalVolume * > useDepth2
A DDSolid represents the shape of a part.
static TrackerG4SimHitNumberingScheme & numberingScheme(const DDCompactView &cpv, const GeometricDet &det)
virtual uint16_t getDepth(G4Step *)
void addLevel(const std::string &name, const int ©Number)
EcalBaseNumber theBaseNumber
bool next()
set current node to the next node in the filtered tree
std::vector< std::pair< unsigned int, DDValue > > DDsvalues_type
std::maps an index to a DDValue. The index corresponds to the index assigned to the name of the std::...
double crystalLength(G4LogicalVolume *)
virtual uint16_t getLayerIDForTimeSim(G4Step *)
DDSolidShape shape(void) const
The type of the solid.
virtual int getTrackID(G4Track *)
double getAttenuation(G4Step *aStep, double birk1, double birk2, double birk3)
void initMap(G4String, const DDCompactView &)
const std::vector< std::string > & strings() const
a reference to the std::string-valued values stored in the given instance of DDValue ...
double curve_LY(G4Step *)
G4StepPoint * preStepPoint
std::vector< double > getDDDArray(const std::string &, const DDsvalues_type &)
virtual uint16_t getRadiationLength(G4Step *)
double getBirkL3(G4Step *)
double calcLightCollectionEfficiencyWeighted(DetId id, double z)
DDsvalues_type mergedSpecifics() const
ECalSD(G4String, const DDCompactView &, const SensitiveDetectorCatalog &, edm::ParameterSet const &p, const SimTrackManager *)
double getResponseWt(G4Track *)
bool firstChild()
set the current node to the first child ...
virtual uint32_t setDetUnitId(G4Step *)
EcalNumberingScheme * numberingScheme
EnergyResolutionVsLumi ageing
std::map< G4LogicalVolume *, double > xtalLMap
void setNumberingScheme(EcalNumberingScheme *)
void setSize(const int &size)
const std::string & name() const
Returns the name.
const DDMaterial & material(void) const
Returns a reference object of the material this LogicalPart is made of.
void setCriteria(const DDValue &nameVal, DDCompOp, DDLogOp l=DDLogOp::AND, bool asString=true, bool merged=true)
virtual uint32_t getUnitID(const EcalBaseNumber &baseNumber) const =0
The DDGenericFilter is a runtime-parametrized Filter looking on DDSpecifcs.
G4ThreeVector setToLocal(const G4ThreeVector &, const G4VTouchable *)