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)
36 return std::find(v.begin(), v.end(), what) != v.end();
42 CaloSD(name, cpv, clg, p, manager,
43 p.getParameter<edm::
ParameterSet>(
"ECalSD").getParameter<int>(
"TimeSliceUnit"),
44 p.getParameter<edm::
ParameterSet>(
"ECalSD").getParameter<bool>(
"IgnoreTrackID")),
86 std::vector<double> tempD =
getDDDArray(
"EnergyWeight",sv);
87 if (tempD.size() > 0) {
if (tempD[0] < 0.1) useWeight =
false; }
96 if (nullNS) scheme = 0;
99 else if (name ==
"EcalHitsES") {
103 }
else {
edm::LogWarning(
"EcalSim") <<
"ECalSD: ReadoutName not supported\n";}
107 LogDebug(
"EcalSim") <<
"Constructing a ECalSD with name " << GetName();
110 edm::LogInfo(
"EcalSim") <<
"ECalSD:: Use of Birks law is set to "
111 << useBirk <<
" with three constants kB = "
112 << birk1 <<
", C1 = " << birk2 <<
", C2 = "
113 << birk3 <<
"\n Use of L3 parametrization "
114 << useBirkL3 <<
" with slope " << birkSlope
115 <<
" and cut off " << birkCut <<
"\n"
116 <<
" Slope for Light yield is set to "
119 edm::LogInfo(
"EcalSim") <<
"ECalSD:: energy deposit is not corrected "
120 <<
" by Birk or light yield curve";
126 <<
" ions below " <<
kmaxIon <<
" MeV\n"
130 if (useWeight)
initMap(name,cpv);
144 G4Track*
theTrack = aStep->GetTrack();
145 double wt2 = theTrack->GetWeight();
146 G4String nameVolume =
preStepPoint->GetPhysicalVolume()->GetName();
153 int pdg = theTrack->GetDefinition()->GetPDGEncoding();
155 double ke = theTrack->GetKineticEnergy()/
MeV;
156 if (((pdg/1000000000 == 1 && ((pdg/10000)%100) > 0 &&
157 ((pdg/10)%100) > 0)) && (ke<
kmaxIon)) weight = 0;
158 if ((pdg == 2212) && (ke <
kmaxProton)) weight = 0;
159 if ((pdg == 2112) && (ke <
kmaxNeutron)) weight = 0;
162 LogDebug(
"EcalSim") <<
"Ignore Track " << theTrack->GetTrackID()
163 <<
" Type " << theTrack->GetDefinition()->GetParticleName()
164 <<
" Kinetic Energy " << ke <<
" MeV";
169 G4LogicalVolume* lv = aStep->GetPreStepPoint()->GetTouchable()->GetVolume(0)->GetLogicalVolume();
178 double edep = aStep->GetTotalEnergyDeposit()*weight*wt1;
195 if(wt2 > 0.0) { edep *= wt2; }
197 LogDebug(
"EcalSim") <<
"ECalSD:: " << nameVolume
198 <<
" Light Collection Efficiency " <<weight <<
":" <<wt1
199 <<
" Weighted Energy Deposit " << edep/
MeV <<
" MeV";
210 G4LogicalVolume* lv =
preStepPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
219 primaryID = aTrack->GetTrackID();
227 G4LogicalVolume* lv = aStep->GetPreStepPoint()->GetTouchable()->GetVolume(0)->GetLogicalVolume();
233 LogDebug(
"EcalSim") <<
"Volume " << lv->GetName() <<
" Depth " <<
ret;
242 G4StepPoint* hitPoint = aStep->GetPreStepPoint();
243 G4LogicalVolume* lv = hitPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
246 G4ThreeVector localPoint =
setToLocal(hitPoint->GetPosition(),
247 hitPoint->GetTouchable());
249 double radl = hitPoint->GetMaterial()->GetRadlen();
250 double detz = (float)(0.5*crlength + localPoint.z());
251 thisX0 = (uint16_t)floor(detz/radl);
268 edm::LogInfo(
"EcalSim") <<
"EcalSD: updates numbering scheme for "
269 << GetName() <<
"\n";
278 G4String attribute =
"ReadOutName";
286 std::vector<G4LogicalVolume*> lvused;
287 const G4LogicalVolumeStore * lvs = G4LogicalVolumeStore::GetInstance();
288 std::vector<G4LogicalVolume *>::const_iterator lvcite;
289 std::map<std::string, G4LogicalVolume *> nameMap;
290 for (
auto lvi = lvs->begin(), lve = lvs->end(); lvi != lve; ++lvi)
291 nameMap.insert(std::make_pair((*lvi)->GetName(), *lvi));
297 G4LogicalVolume* lv = nameMap[lvname];
299 if (strncmp(lvname.c_str(),
depth1Name.c_str(), 4) == 0) {
303 LogDebug(
"EcalSim") <<
"ECalSD::initMap Logical Volume " << lvname
304 <<
" in Depth 1 volume list";
307 G4LogicalVolume* lvr = nameMap[lvname +
"_refl"];
311 LogDebug(
"EcalSim") <<
"ECalSD::initMap Logical Volume " << lvname <<
"_refl"
312 <<
" in Depth 1 volume list";
318 if (strncmp(lvname.c_str(),
depth2Name.c_str(), 4) == 0) {
322 LogDebug(
"EcalSim") <<
"ECalSD::initMap Logical Volume " << lvname
323 <<
" in Depth 2 volume list";
326 G4LogicalVolume* lvr = nameMap[lvname +
"_refl"];
330 LogDebug(
"EcalSim") <<
"ECalSD::initMap Logical Volume " << lvname <<
"_refl"
331 <<
" in Depth 2 volume list";
338 if (!
any(lvused,lv)) {
339 lvused.push_back(lv);
341 const std::vector<double> & paras = sol.
parameters();
343 LogDebug(
"EcalSim") <<
"ECalSD::initMap (for " << sd <<
"): Solid "
344 << lvname <<
" Shape " << sol.
shape()
345 <<
" Parameter 0 = "<< paras[0]
346 <<
" Logical Volume " << lv;
349 double dz = 2*paras[0];
350 xtalLMap.insert(std::pair<G4LogicalVolume*,double>(lv,dz));
351 lv = nameMap[lvname +
"_refl"];
353 xtalLMap.insert(std::pair<G4LogicalVolume*,double>(lv,dz));
360 LogDebug(
"EcalSim") <<
"ECalSD::initMap Logical Volume " << lvname
361 <<
" Material " << matname <<
" in noWeight list";
364 lv = nameMap[lvname];
368 LogDebug(
"EcalSim") <<
"ECalSD::initMap Logical Volume " << lvname
369 <<
" Material " << matname <<
" in noWeight list";
377 LogDebug(
"EcalSim") <<
"ECalSD: Length Table for " << attribute <<
" = "
379 std::map<G4LogicalVolume*,double>::const_iterator ite =
xtalLMap.begin();
381 for (; ite !=
xtalLMap.end(); ite++, i++) {
382 G4String
name =
"Unknown";
383 if (ite->first != 0) name = (ite->first)->GetName();
384 LogDebug(
"EcalSim") <<
" " << i <<
" " << ite->first <<
" " << name
385 <<
" L = " << ite->second;
392 G4StepPoint* stepPoint = aStep->GetPreStepPoint();
393 G4LogicalVolume* lv = stepPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
396 G4ThreeVector localPoint =
setToLocal(stepPoint->GetPosition(),
397 stepPoint->GetTouchable());
403 double depth = 0.5 * crlength + localPoint.z();
405 if (depth >= -0.1 || depth <= crlength+0.1)
409 double dapd = 0.5 * crlength - localPoint.z();
410 if (dapd >= -0.1 || dapd <= crlength+0.1) {
414 edm::LogWarning(
"EcalSim") <<
"ECalSD: light coll curve : wrong distance "
415 <<
"to APD " << dapd <<
" crlength = "
416 << crlength <<
" crystal name = " <<lv->GetName()
417 <<
" z of localPoint = " << localPoint.z()
418 <<
" take weight = " <<
weight;
422 LogDebug(
"EcalSim") <<
"ECalSD, light coll curve : " << dapd
423 <<
" crlength = " << crlength
424 <<
" crystal name = " << lv->GetName()
425 <<
" z of localPoint = " << localPoint.z()
426 <<
" take weight = " <<
weight;
434 std::map<G4LogicalVolume*,double>::const_iterator ite =
xtalLMap.find(lv);
435 if (ite !=
xtalLMap.end()) length = ite->second;
442 const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
443 int theSize = touch->GetHistoryDepth()+1;
447 for (
int ii = 0;
ii < theSize ;
ii++) {
450 LogDebug(
"EcalSim") <<
"ECalSD::getBaseNumber(): Adding level " <<
ii
451 <<
": " << touch->GetVolume(
ii)->GetName() <<
"["
452 << touch->GetReplicaNumber(
ii) <<
"]";
462 double charge = aStep->GetPreStepPoint()->GetCharge();
464 if (charge != 0. && aStep->GetStepLength() > 0) {
465 G4Material* mat = aStep->GetPreStepPoint()->GetMaterial();
466 double density = mat->GetDensity();
467 double dedx = aStep->GetTotalEnergyDeposit()/aStep->GetStepLength();
468 double rkb =
birk1/density;
472 else if (weight > 1.) weight = 1.;
475 LogDebug(
"EcalSim") <<
"ECalSD::getBirkL3 in " << mat->GetName()
476 <<
" Charge " << charge <<
" dE/dx " << dedx
477 <<
" Birk Const " << rkb <<
" Weight = " << weight
478 <<
" dE " << aStep->GetTotalEnergyDeposit();
489 LogDebug(
"EcalSim") <<
"ECalSD:getDDDArray called for " << str;
496 const std::vector<double> & fvec = value.
doubles();
499 std::vector<double> fvec;
508 LogDebug(
"EcalSim") <<
"ECalSD:getStringArray called for " << str;
515 const std::vector<std::string> & fvec = value.
strings();
518 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
void addFilter(const DDFilter &, log_op op=AND)
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 *)
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 *)
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 ...
void setCriteria(const DDValue &nameVal, comp_op, log_op l=AND, bool asString=true, bool merged=true)
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.
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 *)