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")),
84 std::vector<double> tempD =
getDDDArray(
"EnergyWeight",sv);
85 if (tempD.size() > 0) {
if (tempD[0] < 0.1) useWeight =
false; }
94 if (nullNS) scheme = 0;
97 else if (name ==
"EcalHitsES") {
101 }
else {
edm::LogWarning(
"EcalSim") <<
"ECalSD: ReadoutName not supported\n";}
105 LogDebug(
"EcalSim") <<
"Constructing a ECalSD with name " << GetName();
108 edm::LogInfo(
"EcalSim") <<
"ECalSD:: Use of Birks law is set to "
109 << useBirk <<
" with three constants kB = "
110 << birk1 <<
", C1 = " << birk2 <<
", C2 = "
111 << birk3 <<
"\n Use of L3 parametrization "
112 << useBirkL3 <<
" with slope " << birkSlope
113 <<
" and cut off " << birkCut <<
"\n"
114 <<
" Slope for Light yield is set to "
117 edm::LogInfo(
"EcalSim") <<
"ECalSD:: energy deposit is not corrected "
118 <<
" by Birk or light yield curve";
124 <<
" ions below " <<
kmaxIon <<
" MeV\n"
128 if (useWeight)
initMap(name,cpv);
142 G4Track*
theTrack = aStep->GetTrack();
143 double wt2 = theTrack->GetWeight();
144 G4String nameVolume =
preStepPoint->GetPhysicalVolume()->GetName();
151 int pdg = theTrack->GetDefinition()->GetPDGEncoding();
153 double ke = theTrack->GetKineticEnergy()/MeV;
154 if (((pdg/1000000000 == 1 && ((pdg/10000)%100) > 0 &&
155 ((pdg/10)%100) > 0)) && (ke<
kmaxIon)) weight = 0;
156 if ((pdg == 2212) && (ke <
kmaxProton)) weight = 0;
157 if ((pdg == 2112) && (ke <
kmaxNeutron)) weight = 0;
160 LogDebug(
"EcalSim") <<
"Ignore Track " << theTrack->GetTrackID()
161 <<
" Type " << theTrack->GetDefinition()->GetParticleName()
162 <<
" Kinetic Energy " << ke <<
" MeV";
167 G4LogicalVolume* lv = aStep->GetPreStepPoint()->GetTouchable()->GetVolume(0)->GetLogicalVolume();
176 double edep = aStep->GetTotalEnergyDeposit()*weight*wt1;
193 if(wt2 > 0.0) { edep *= wt2; }
195 LogDebug(
"EcalSim") <<
"ECalSD:: " << nameVolume
196 <<
" Light Collection Efficiency " <<weight <<
":" <<wt1
197 <<
" Weighted Energy Deposit " << edep/MeV <<
" MeV";
208 G4LogicalVolume* lv =
preStepPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
217 primaryID = aTrack->GetTrackID();
225 G4LogicalVolume* lv = aStep->GetPreStepPoint()->GetTouchable()->GetVolume(0)->GetLogicalVolume();
231 LogDebug(
"EcalSim") <<
"Volume " << lv->GetName() <<
" Depth " <<
ret;
240 G4StepPoint* hitPoint = aStep->GetPreStepPoint();
241 G4LogicalVolume* lv = hitPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
244 G4ThreeVector localPoint =
setToLocal(hitPoint->GetPosition(),
245 hitPoint->GetTouchable());
247 double radl = hitPoint->GetMaterial()->GetRadlen();
248 double detz = (float)(0.5*crlength + localPoint.z());
249 thisX0 = (uint16_t)floor(detz/radl);
266 edm::LogInfo(
"EcalSim") <<
"EcalSD: updates numbering scheme for "
267 << GetName() <<
"\n";
276 G4String attribute =
"ReadOutName";
284 std::vector<G4LogicalVolume*> lvused;
285 const G4LogicalVolumeStore * lvs = G4LogicalVolumeStore::GetInstance();
286 std::vector<G4LogicalVolume *>::const_iterator lvcite;
287 std::map<std::string, G4LogicalVolume *> nameMap;
288 for (
auto lvi = lvs->begin(), lve = lvs->end(); lvi != lve; ++lvi)
289 nameMap.insert(std::make_pair((*lvi)->GetName(), *lvi));
295 G4LogicalVolume* lv = nameMap[lvname];
297 if (strncmp(lvname.c_str(),
depth1Name.c_str(), 4) == 0) {
301 LogDebug(
"EcalSim") <<
"ECalSD::initMap Logical Volume " << lvname
302 <<
" in Depth 1 volume list";
305 G4LogicalVolume* lvr = nameMap[lvname +
"_refl"];
309 LogDebug(
"EcalSim") <<
"ECalSD::initMap Logical Volume " << lvname <<
"_refl"
310 <<
" in Depth 1 volume list";
316 if (strncmp(lvname.c_str(),
depth2Name.c_str(), 4) == 0) {
320 LogDebug(
"EcalSim") <<
"ECalSD::initMap Logical Volume " << lvname
321 <<
" in Depth 2 volume list";
324 G4LogicalVolume* lvr = nameMap[lvname +
"_refl"];
328 LogDebug(
"EcalSim") <<
"ECalSD::initMap Logical Volume " << lvname <<
"_refl"
329 <<
" in Depth 2 volume list";
336 if (!
any(lvused,lv)) {
337 lvused.push_back(lv);
339 const std::vector<double> & paras = sol.
parameters();
341 LogDebug(
"EcalSim") <<
"ECalSD::initMap (for " << sd <<
"): Solid "
342 << lvname <<
" Shape " << sol.
shape()
343 <<
" Parameter 0 = "<< paras[0]
344 <<
" Logical Volume " << lv;
347 double dz = 2*paras[0];
348 xtalLMap.insert(std::pair<G4LogicalVolume*,double>(lv,dz));
349 lv = nameMap[lvname +
"_refl"];
351 xtalLMap.insert(std::pair<G4LogicalVolume*,double>(lv,dz));
358 LogDebug(
"EcalSim") <<
"ECalSD::initMap Logical Volume " << lvname
359 <<
" Material " << matname <<
" in noWeight list";
362 lv = nameMap[lvname];
366 LogDebug(
"EcalSim") <<
"ECalSD::initMap Logical Volume " << lvname
367 <<
" Material " << matname <<
" in noWeight list";
375 LogDebug(
"EcalSim") <<
"ECalSD: Length Table for " << attribute <<
" = "
377 std::map<G4LogicalVolume*,double>::const_iterator ite =
xtalLMap.begin();
379 for (; ite !=
xtalLMap.end(); ite++, i++) {
380 G4String
name =
"Unknown";
381 if (ite->first != 0) name = (ite->first)->GetName();
382 LogDebug(
"EcalSim") <<
" " << i <<
" " << ite->first <<
" " << name
383 <<
" L = " << ite->second;
390 G4StepPoint* stepPoint = aStep->GetPreStepPoint();
391 G4LogicalVolume* lv = stepPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
394 G4ThreeVector localPoint =
setToLocal(stepPoint->GetPosition(),
395 stepPoint->GetTouchable());
401 double depth = 0.5 * crlength + localPoint.z();
403 if (depth >= -0.1 || depth <= crlength+0.1)
407 double dapd = 0.5 * crlength - localPoint.z();
408 if (dapd >= -0.1 || dapd <= crlength+0.1) {
412 edm::LogWarning(
"EcalSim") <<
"ECalSD: light coll curve : wrong distance "
413 <<
"to APD " << dapd <<
" crlength = "
414 << crlength <<
" crystal name = " <<lv->GetName()
415 <<
" z of localPoint = " << localPoint.z()
416 <<
" take weight = " <<
weight;
420 LogDebug(
"EcalSim") <<
"ECalSD, light coll curve : " << dapd
421 <<
" crlength = " << crlength
422 <<
" crystal name = " << lv->GetName()
423 <<
" z of localPoint = " << localPoint.z()
424 <<
" take weight = " <<
weight;
432 std::map<G4LogicalVolume*,double>::const_iterator ite =
xtalLMap.find(lv);
433 if (ite !=
xtalLMap.end()) length = ite->second;
440 const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
441 int theSize = touch->GetHistoryDepth()+1;
445 for (
int ii = 0;
ii < theSize ;
ii++) {
448 LogDebug(
"EcalSim") <<
"ECalSD::getBaseNumber(): Adding level " <<
ii
449 <<
": " << touch->GetVolume(
ii)->GetName() <<
"["
450 << touch->GetReplicaNumber(
ii) <<
"]";
460 double charge = aStep->GetPreStepPoint()->GetCharge();
462 if (charge != 0. && aStep->GetStepLength() > 0) {
463 G4Material* mat = aStep->GetPreStepPoint()->GetMaterial();
464 double density = mat->GetDensity();
465 double dedx = aStep->GetTotalEnergyDeposit()/aStep->GetStepLength();
466 double rkb =
birk1/density;
470 else if (weight > 1.) weight = 1.;
473 LogDebug(
"EcalSim") <<
"ECalSD::getBirkL3 in " << mat->GetName()
474 <<
" Charge " << charge <<
" dE/dx " << dedx
475 <<
" Birk Const " << rkb <<
" Weight = " << weight
476 <<
" dE " << aStep->GetTotalEnergyDeposit();
487 LogDebug(
"EcalSim") <<
"ECalSD:getDDDArray called for " << str;
494 const std::vector<double> & fvec = value.
doubles();
497 std::vector<double> fvec;
506 LogDebug(
"EcalSim") <<
"ECalSD:getStringArray called for " << str;
513 const std::vector<std::string> & fvec = value.
strings();
516 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 setLumies(double x, double y)
void getBaseNumber(const G4Step *)
virtual double getEnergyDeposit(G4Step *)
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
double getResponseWt(G4Track *)
ECalSD(G4String, const DDCompactView &, SensitiveDetectorCatalog &, edm::ParameterSet const &p, const SimTrackManager *)
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.