24 #include "G4LogicalVolumeStore.hh" 25 #include "G4LogicalVolume.hh" 28 #include "G4VProcess.hh" 30 #include "G4SystemOfUnits.hh" 37 bool any(
const std::vector<T> &
v,
const T &what) {
38 return std::find(v.begin(), v.end(), what) != v.end();
44 CaloSD(name, cpv, clg, p, manager,
90 if (!tempD.empty()) {
if (tempD[0] < 0.1) useWeight =
false; }
92 edm::LogInfo(
"EcalSim") <<
"ECalSD:: useWeight " << tempD.size() <<
":" 93 << useWeight << std::endl;
108 }
else if (name ==
"EcalHitsEB") {
111 }
else if (name ==
"EcalHitsEE") {
114 }
else if (name ==
"EcalHitsES") {
124 edm::LogInfo(
"EcalSim") <<
"Constructing a ECalSD with name " << GetName();
127 edm::LogInfo(
"EcalSim") <<
"ECalSD:: Use of Birks law is set to " 128 << useBirk <<
" with three constants kB = " 129 << birk1 <<
", C1 = " << birk2 <<
", C2 = " 130 << birk3 <<
"\n Use of L3 parametrization " 131 << useBirkL3 <<
" with slope " << birkSlope
132 <<
" and cut off " << birkCut <<
"\n" 133 <<
" Slope for Light yield is set to " 136 edm::LogInfo(
"EcalSim") <<
"ECalSD:: energy deposit is not corrected " 137 <<
" by Birk or light yield curve";
141 <<
"\tprotons below " <<
kmaxProton <<
" MeV," 143 <<
"\tions below " <<
kmaxIon <<
" MeV" 146 <<
"\n\tstoreRL" << storeRL <<
":" << scaleRL
149 if (useWeight)
initMap(name,cpv);
154 static const std::string ctype[4] = {
"EB",
"EBref",
"EE",
"EERef"};
155 for (
int k=0;
k<4; ++
k) {
158 double xmin = (
k > 1) ? 3000.0 : 1000.0;
159 g2L_[
k] = ecDir.
make<TH2F>(name.c_str(),title.c_str(),100,
xmin,
160 xmin+1000.,100,0.0,3000.);
163 for (
int k=0;
k<4; ++
k) g2L_[
k] = 0;
175 if (aStep ==
nullptr) {
179 G4Track*
theTrack = aStep->GetTrack();
180 double wt2 = theTrack->GetWeight();
181 G4String nameVolume =
preStepPoint->GetPhysicalVolume()->GetName();
188 int pdg = theTrack->GetDefinition()->GetPDGEncoding();
190 double ke = theTrack->GetKineticEnergy()/
MeV;
191 if (((pdg/1000000000 == 1 && ((pdg/10000)%100) > 0 &&
192 ((pdg/10)%100) > 0)) && (ke<
kmaxIon)) weight = 0;
193 if ((pdg == 2212) && (ke <
kmaxProton)) weight = 0;
194 if ((pdg == 2112) && (ke <
kmaxNeutron)) weight = 0;
197 edm::LogInfo(
"EcalSim") <<
"Ignore Track " << theTrack->GetTrackID()
198 <<
" Type " << theTrack->GetDefinition()->GetParticleName()
199 <<
" Kinetic Energy " << ke <<
" MeV";
204 G4LogicalVolume* lv = aStep->GetPreStepPoint()->GetTouchable()->GetVolume(0)->GetLogicalVolume();
213 double edep = aStep->GetTotalEnergyDeposit()*weight*wt1;
235 if(wt2 > 0.0) { edep *= wt2; }
238 <<
" Light Collection Efficiency " << weight <<
":" 239 << wt1 <<
" Weighted Energy Deposit " << edep/
MeV 251 G4LogicalVolume* lv =
preStepPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
260 primaryID = aTrack->GetTrackID();
268 G4LogicalVolume* lv = aStep->GetPreStepPoint()->GetTouchable()->GetVolume(0)->GetLogicalVolume();
272 uint16_t depth1 = (ite ==
xtalLMap.end()) ? 0 : (((ite->second) >= 0) ? 0 :
281 edm::LogVerbatim(
"EcalSim") <<
"ECalSD::Depth " << std::hex << depth1 <<
":" 282 << depth2 <<
":" << depth <<
std::dec <<
" L " 283 << (ite ==
xtalLMap.end()) <<
":" <<ite->second;
291 if (aStep !=
nullptr) {
292 G4StepPoint* hitPoint = aStep->GetPreStepPoint();
293 G4LogicalVolume* lv = hitPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
298 G4ThreeVector localPoint =
setToLocal(hitPoint->GetPosition(),
299 hitPoint->GetTouchable());
300 double radl = hitPoint->GetMaterial()->GetRadlen();
302 thisX0 = (uint16_t)floor(
scaleRL*depth/radl);
305 int k1 = (lvname.find(
"EFRY")!=std::string::npos) ? 2 : 0;
306 int k2 = (lvname.find(
"refl")!=std::string::npos) ? 1 : 0;
308 double rz = (k1 == 0) ? (hitPoint->GetPosition()).
rho() :
311 << kk <<
" rz " << rz <<
" D " << thisX0;
312 g2L_[
kk]->Fill(rz,thisX0);
317 << hitPoint->GetPosition() <<
":" 318 << (hitPoint->GetPosition()).
rho()
319 <<
" Local " << localPoint
320 <<
" Crystal Length " << crlength
321 <<
" Radl " << radl <<
" DetZ " << detz
322 <<
" Index " << thisX0
332 float layerSize = 1*cm;
335 if (aStep !=
nullptr ) {
336 G4StepPoint* hitPoint = aStep->GetPreStepPoint();
337 G4LogicalVolume* lv = hitPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
338 G4ThreeVector localPoint =
setToLocal(hitPoint->GetPosition(),
339 hitPoint->GetTouchable());
342 return (
int)detz/layerSize;
357 if (scheme !=
nullptr) {
358 edm::LogInfo(
"EcalSim") <<
"EcalSD: updates numbering scheme for " 368 G4String attribute =
"ReadOutName";
373 std::vector<G4LogicalVolume*> lvused;
374 const G4LogicalVolumeStore * lvs = G4LogicalVolumeStore::GetInstance();
375 std::vector<G4LogicalVolume *>::const_iterator lvcite;
376 std::map<std::string, G4LogicalVolume *> nameMap;
377 for (
auto lvi = lvs->begin(), lve = lvs->end(); lvi != lve; ++lvi)
378 nameMap.insert(std::make_pair((*lvi)->GetName(), *lvi));
384 G4LogicalVolume* lv = nameMap[lvname];
385 int ibec = (lvname.find(
"EFRY") == std::string::npos) ? 0 : 1;
386 int iref = (lvname.find(
"refl") == std::string::npos) ? 0 : 1;
387 int type = (ibec+iref == 1) ? 1 : -1;
389 if (strncmp(lvname.c_str(),
depth1Name.c_str(), 4) == 0) {
393 edm::LogInfo(
"EcalSim") <<
"ECalSD::initMap Logical Volume " 394 << lvname <<
" in Depth 1 volume list";
397 G4LogicalVolume* lvr = nameMap[lvname +
"_refl"];
401 edm::LogInfo(
"EcalSim") <<
"ECalSD::initMap Logical Volume " 403 <<
" in Depth 1 volume list";
409 if (strncmp(lvname.c_str(),
depth2Name.c_str(), 4) == 0) {
413 edm::LogInfo(
"EcalSim") <<
"ECalSD::initMap Logical Volume " 414 << lvname <<
" in Depth 2 volume list";
417 G4LogicalVolume* lvr = nameMap[lvname +
"_refl"];
421 edm::LogInfo(
"EcalSim") <<
"ECalSD::initMap Logical Volume " 423 <<
" in Depth 2 volume list";
430 if (!
any(lvused,lv)) {
431 lvused.push_back(lv);
433 const std::vector<double> & paras = sol.
parameters();
435 edm::LogInfo(
"EcalSim") <<
"ECalSD::initMap (for " << sd
436 <<
"): Solid " << lvname <<
" Shape " 437 << sol.
shape() <<
" Parameter 0 = " 438 << paras[0] <<
" Logical Volume " << lv;
441 double dz = 2*paras[0];
442 xtalLMap.insert(std::pair<G4LogicalVolume*,double>(lv,dz*type));
443 lv = nameMap[lvname +
"_refl"];
445 xtalLMap.insert(std::pair<G4LogicalVolume*,double>(lv,-dz*type));
453 edm::LogInfo(
"EcalSim") <<
"ECalSD::initMap Logical Volume " 454 << lvname <<
" Material " << matname
455 <<
" in noWeight list";
458 lv = nameMap[lvname];
462 edm::LogInfo(
"EcalSim") <<
"ECalSD::initMap Logical Volume " 463 << lvname <<
" Material " << matname
464 <<
" in noWeight list";
472 edm::LogInfo(
"EcalSim") <<
"ECalSD: Length Table for " << attribute <<
" = " 476 G4String
name(
"Unknown");
477 if (ite.first != 0) name = (ite.first)->GetName();
478 edm::LogInfo(
"EcalSim") <<
" " << i <<
" " << ite.first <<
" " << name
479 <<
" L = " << ite.second;
487 G4StepPoint* stepPoint = aStep->GetPreStepPoint();
488 G4LogicalVolume* lv = stepPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
491 G4ThreeVector localPoint =
setToLocal(stepPoint->GetPosition(),
492 stepPoint->GetTouchable());
499 if (depth >= -0.1 || depth <= crlength+0.1)
502 double dapd = crlength -
depth;
503 if (dapd >= -0.1 || dapd <= crlength+0.1) {
507 edm::LogWarning(
"EcalSim") <<
"ECalSD: light coll curve : wrong distance " 508 <<
"to APD " << dapd <<
" crlength = " 509 << crlength <<
" crystal name = " <<lv->GetName()
510 <<
" z of localPoint = " << localPoint.z()
511 <<
" take weight = " <<
weight;
525 const G4ThreeVector& localPoint) {
529 (
std::abs(0.5*(ite->second)+localPoint.z()));
537 const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
538 int theSize = touch->GetHistoryDepth()+1;
542 for (
int ii = 0;
ii < theSize ;
ii++) {
545 edm::LogInfo(
"EcalSim") <<
"ECalSD::getBaseNumber(): Adding level " <<
ii 546 <<
": " << touch->GetVolume(
ii)->GetName() <<
"[" 547 << touch->GetReplicaNumber(
ii) <<
"]";
557 double charge = aStep->GetPreStepPoint()->GetCharge();
559 if (charge != 0. && aStep->GetStepLength() > 0) {
560 G4Material* mat = aStep->GetPreStepPoint()->GetMaterial();
561 double density = mat->GetDensity();
562 double dedx = aStep->GetTotalEnergyDeposit()/aStep->GetStepLength();
563 double rkb =
birk1/density;
567 else if (weight > 1.) weight = 1.;
570 edm::LogInfo(
"EcalSim") <<
"ECalSD::getBirkL3 in " << mat->GetName()
571 <<
" Charge " << charge <<
" dE/dx " << dedx
572 <<
" Birk Const " << rkb <<
" Weight = " << weight
573 <<
" dE " << aStep->GetTotalEnergyDeposit();
591 const std::vector<double> & fvec = value.
doubles();
594 std::vector<double> fvec;
610 const std::vector<std::string> & fvec = value.
strings();
613 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 *)
static const int kEcalDepthRefz
bool any(const std::vector< T > &v, const T &what)
uint16_t getDepth(G4Step *) override
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)
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::...
static const int kEcalDepthMask
double crystalLength(G4LogicalVolume *)
Abs< T >::type abs(const T &t)
virtual uint16_t getLayerIDForTimeSim(G4Step *)
DDSolidShape shape(void) const
The type of the solid.
static const int kEcalDepthOffset
T * make(const Args &...args) const
make new ROOT object
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
int getTrackID(G4Track *) override
TFileDirectory mkdir(const std::string &dir, const std::string &descr="")
create a new subdirectory
virtual uint32_t getUnitID(const EcalBaseNumber &baseNumber) const =0
double crystalDepth(G4LogicalVolume *, const G4ThreeVector &)
std::vector< double > getDDDArray(const std::string &, const DDsvalues_type &)
virtual uint16_t getRadiationLength(G4Step *)
uint32_t setDetUnitId(G4Step *) override
double getBirkL3(G4Step *)
double calcLightCollectionEfficiencyWeighted(DetId id, double z)
DDsvalues_type mergedSpecifics() const
double getEnergyDeposit(G4Step *) override
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 ...
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.
G4ThreeVector setToLocal(const G4ThreeVector &, const G4VTouchable *)