21 #include "G4LogicalVolumeStore.hh"
22 #include "G4LogicalVolume.hh"
25 #include "G4VProcess.hh"
32 bool any(
const std::vector<T> &
v,
const T &what)
34 return std::find(v.begin(), v.end(), what) != v.end();
40 CaloSD(name, cpv, clg, p, manager,
41 p.getParameter<edm::
ParameterSet>(
"ECalSD").getParameter<int>(
"TimeSliceUnit"),
42 p.getParameter<edm::
ParameterSet>(
"ECalSD").getParameter<bool>(
"IgnoreTrackID")),
70 std::string attribute =
"ReadOutName";
80 std::vector<double> tempD =
getDDDArray(
"EnergyWeight",sv);
81 if (tempD.size() > 0) {
if (tempD[0] < 0.1) useWeight =
false; }
90 if (nullNS) scheme = 0;
93 else if (name ==
"EcalHitsES") {
97 }
else {
edm::LogWarning(
"EcalSim") <<
"ECalSD: ReadoutName not supported\n";}
101 LogDebug(
"EcalSim") <<
"Constructing a ECalSD with name " << GetName();
104 edm::LogInfo(
"EcalSim") <<
"ECalSD:: Use of Birks law is set to "
105 << useBirk <<
" with three constants kB = "
106 << birk1 <<
", C1 = " << birk2 <<
", C2 = "
107 << birk3 <<
"\n Use of L3 parametrization "
108 << useBirkL3 <<
" with slope " << birkSlope
109 <<
" and cut off " << birkCut <<
"\n"
110 <<
" Slope for Light yield is set to "
113 edm::LogInfo(
"EcalSim") <<
"ECalSD:: energy deposit is not corrected "
114 <<
" by Birk or light yield curve";
120 <<
" ions below " <<
kmaxIon <<
" MeV\n"
124 if (useWeight)
initMap(name,cpv);
138 G4String nameVolume =
preStepPoint->GetPhysicalVolume()->GetName();
143 G4Track*
theTrack = aStep->GetTrack();
146 int pdg = theTrack->GetDefinition()->GetPDGEncoding();
148 double ke = theTrack->GetKineticEnergy()/MeV;
149 if (((pdg/1000000000 == 1 && ((pdg/10000)%100) > 0 &&
150 ((pdg/10)%100) > 0)) && (ke<
kmaxIon)) weight = 0;
151 if ((pdg == 2212) && (ke <
kmaxProton)) weight = 0;
152 if ((pdg == 2112) && (ke <
kmaxNeutron)) weight = 0;
155 LogDebug(
"EcalSim") <<
"Ignore Track " << theTrack->GetTrackID()
156 <<
" Type " << theTrack->GetDefinition()->GetParticleName()
157 <<
" Kinetic Energy " << ke <<
" MeV";
162 G4LogicalVolume* lv = aStep->GetPreStepPoint()->GetTouchable()->GetVolume(0)->GetLogicalVolume();
171 double edep = aStep->GetTotalEnergyDeposit() * weight * wt1;
173 LogDebug(
"EcalSim") <<
"ECalSD:: " << nameVolume
174 <<
" Light Collection Efficiency " <<weight <<
":" <<wt1
175 <<
" Weighted Energy Deposit " << edep/MeV <<
" MeV";
186 G4LogicalVolume* lv =
preStepPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
195 primaryID = aTrack->GetTrackID();
203 G4LogicalVolume* lv = aStep->GetPreStepPoint()->GetTouchable()->GetVolume(0)->GetLogicalVolume();
209 LogDebug(
"EcalSim") <<
"Volume " << lv->GetName() <<
" Depth " <<
ret;
218 G4StepPoint* hitPoint = aStep->GetPreStepPoint();
219 G4LogicalVolume* lv = hitPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
222 G4ThreeVector localPoint =
setToLocal(hitPoint->GetPosition(),
223 hitPoint->GetTouchable());
225 double radl = hitPoint->GetMaterial()->GetRadlen();
226 double detz = (float)(0.5*crlength + localPoint.z());
227 thisX0 = (uint16_t)floor(detz/radl);
244 edm::LogInfo(
"EcalSim") <<
"EcalSD: updates numbering scheme for "
245 << GetName() <<
"\n";
254 G4String attribute =
"ReadOutName";
262 std::vector<G4LogicalVolume*> lvused;
263 const G4LogicalVolumeStore * lvs = G4LogicalVolumeStore::GetInstance();
264 std::vector<G4LogicalVolume *>::const_iterator lvcite;
265 std::map<std::string, G4LogicalVolume *> nameMap;
266 for (
auto lvi = lvs->begin(), lve = lvs->end(); lvi != lve; ++lvi)
267 nameMap.insert(std::make_pair((*lvi)->GetName(), *lvi));
273 G4LogicalVolume* lv = nameMap[lvname];
275 if (strncmp(lvname.c_str(),
depth1Name.c_str(), 4) == 0) {
279 LogDebug(
"EcalSim") <<
"ECalSD::initMap Logical Volume " << lvname
280 <<
" in Depth 1 volume list";
283 G4LogicalVolume* lvr = nameMap[lvname +
"_refl"];
287 LogDebug(
"EcalSim") <<
"ECalSD::initMap Logical Volume " << lvname <<
"_refl"
288 <<
" in Depth 1 volume list";
294 if (strncmp(lvname.c_str(),
depth2Name.c_str(), 4) == 0) {
298 LogDebug(
"EcalSim") <<
"ECalSD::initMap Logical Volume " << lvname
299 <<
" in Depth 2 volume list";
302 G4LogicalVolume* lvr = nameMap[lvname +
"_refl"];
306 LogDebug(
"EcalSim") <<
"ECalSD::initMap Logical Volume " << lvname <<
"_refl"
307 <<
" in Depth 2 volume list";
314 if (!
any(lvused,lv)) {
315 lvused.push_back(lv);
317 const std::vector<double> & paras = sol.
parameters();
319 LogDebug(
"EcalSim") <<
"ECalSD::initMap (for " << sd <<
"): Solid "
320 << lvname <<
" Shape " << sol.
shape()
321 <<
" Parameter 0 = "<< paras[0]
322 <<
" Logical Volume " << lv;
325 double dz = 2*paras[0];
326 xtalLMap.insert(std::pair<G4LogicalVolume*,double>(lv,dz));
327 lv = nameMap[lvname +
"_refl"];
329 xtalLMap.insert(std::pair<G4LogicalVolume*,double>(lv,dz));
336 LogDebug(
"EcalSim") <<
"ECalSD::initMap Logical Volume " << lvname
337 <<
" Material " << matname <<
" in noWeight list";
340 lv = nameMap[lvname];
344 LogDebug(
"EcalSim") <<
"ECalSD::initMap Logical Volume " << lvname
345 <<
" Material " << matname <<
" in noWeight list";
353 LogDebug(
"EcalSim") <<
"ECalSD: Length Table for " << attribute <<
" = "
355 std::map<G4LogicalVolume*,double>::const_iterator ite =
xtalLMap.begin();
357 for (; ite !=
xtalLMap.end(); ite++, i++) {
358 G4String
name =
"Unknown";
359 if (ite->first != 0) name = (ite->first)->GetName();
360 LogDebug(
"EcalSim") <<
" " << i <<
" " << ite->first <<
" " << name
361 <<
" L = " << ite->second;
368 G4StepPoint* stepPoint = aStep->GetPreStepPoint();
369 G4LogicalVolume* lv = stepPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
372 G4ThreeVector localPoint =
setToLocal(stepPoint->GetPosition(),
373 stepPoint->GetTouchable());
375 double dapd = 0.5 * crlength - localPoint.z();
376 if (dapd >= -0.1 || dapd <= crlength+0.1) {
380 edm::LogWarning(
"EcalSim") <<
"ECalSD: light coll curve : wrong distance "
381 <<
"to APD " << dapd <<
" crlength = "
382 << crlength <<
" crystal name = " <<lv->GetName()
383 <<
" z of localPoint = " << localPoint.z()
384 <<
" take weight = " <<
weight;
387 LogDebug(
"EcalSim") <<
"ECalSD, light coll curve : " << dapd
388 <<
" crlength = " << crlength
389 <<
" crystal name = " << lv->GetName()
390 <<
" z of localPoint = " << localPoint.z()
391 <<
" take weight = " <<
weight;
399 std::map<G4LogicalVolume*,double>::const_iterator ite =
xtalLMap.find(lv);
400 if (ite !=
xtalLMap.end()) length = ite->second;
407 const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
408 int theSize = touch->GetHistoryDepth()+1;
412 for (
int ii = 0; ii < theSize ; ii++) {
415 LogDebug(
"EcalSim") <<
"ECalSD::getBaseNumber(): Adding level " << ii
416 <<
": " << touch->GetVolume(ii)->GetName() <<
"["
417 << touch->GetReplicaNumber(ii) <<
"]";
427 double charge = aStep->GetPreStepPoint()->GetCharge();
429 if (charge != 0. && aStep->GetStepLength() > 0) {
430 G4Material* mat = aStep->GetPreStepPoint()->GetMaterial();
431 double density = mat->GetDensity();
432 double dedx = aStep->GetTotalEnergyDeposit()/aStep->GetStepLength();
433 double rkb =
birk1/density;
437 else if (weight > 1.) weight = 1.;
440 LogDebug(
"EcalSim") <<
"ECalSD::getBirkL3 in " << mat->GetName()
441 <<
" Charge " << charge <<
" dE/dx " << dedx
442 <<
" Birk Const " << rkb <<
" Weight = " << weight
443 <<
" dE " << aStep->GetTotalEnergyDeposit();
454 LogDebug(
"EcalSim") <<
"ECalSD:getDDDArray called for " << str;
461 const std::vector<double> & fvec = value.
doubles();
464 std::vector<double> fvec;
473 LogDebug(
"EcalSim") <<
"ECalSD:getStringArray called for " << str;
480 const std::vector<std::string> & fvec = value.
strings();
483 std::vector<std::string> fvec;
G4ThreeVector setToLocal(G4ThreeVector, const G4VTouchable *)
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 getBaseNumber(const G4Step *)
virtual double getEnergyDeposit(G4Step *)
virtual int getTrackID(G4Track *)
ECalSD(G4String, const DDCompactView &, SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
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 *)
DDsvalues_type mergedSpecifics() const
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
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.