CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
ECalSD Class Reference

#include <ECalSD.h>

Inheritance diagram for ECalSD:
CaloSD SensitiveCaloDetector Observer< const BeginOfRun * > Observer< const BeginOfEvent * > Observer< const BeginOfTrack * > Observer< const EndOfTrack * > Observer< const EndOfEvent * > SensitiveDetector

Public Member Functions

 ECalSD (G4String, const DDCompactView &, SensitiveDetectorCatalog &, edm::ParameterSet const &p, const SimTrackManager *)
 
virtual uint16_t getDepth (G4Step *)
 
virtual double getEnergyDeposit (G4Step *)
 
virtual uint16_t getRadiationLength (G4Step *)
 
virtual int getTrackID (G4Track *)
 
virtual uint32_t setDetUnitId (G4Step *)
 
void setNumberingScheme (EcalNumberingScheme *)
 
virtual ~ECalSD ()
 
- Public Member Functions inherited from CaloSD
 CaloSD (G4String aSDname, const DDCompactView &cpv, SensitiveDetectorCatalog &clg, edm::ParameterSet const &p, const SimTrackManager *, int tSlice=1, bool ignoreTkID=false)
 
virtual void clear ()
 
virtual void DrawAll ()
 
virtual void EndOfEvent (G4HCofThisEvent *eventHC)
 
void fillHits (edm::PCaloHitContainer &, std::string n)
 
virtual void Initialize (G4HCofThisEvent *HCE)
 
virtual void PrintAll ()
 
virtual bool ProcessHits (G4Step *step, G4TouchableHistory *tHistory)
 
virtual bool ProcessHits (G4GFlashSpot *aSpot, G4TouchableHistory *)
 
virtual ~CaloSD ()
 
- Public Member Functions inherited from SensitiveCaloDetector
 SensitiveCaloDetector (std::string &iname, const DDCompactView &cpv, SensitiveDetectorCatalog &clg, edm::ParameterSet const &p)
 
- Public Member Functions inherited from SensitiveDetector
virtual void AssignSD (std::string &vname)
 
Local3DPoint ConvertToLocal3DPoint (const G4ThreeVector &point)
 
Local3DPoint FinalStepPosition (G4Step *s, coordinates)
 
virtual std::vector< std::string > getNames ()
 
Local3DPoint InitialStepPosition (G4Step *s, coordinates)
 
std::string nameOfSD ()
 
void NaNTrap (G4Step *step)
 
void Register ()
 
 SensitiveDetector (std::string &iname, const DDCompactView &cpv, SensitiveDetectorCatalog &, edm::ParameterSet const &p)
 
virtual ~SensitiveDetector ()
 
- Public Member Functions inherited from Observer< const BeginOfRun * >
 Observer ()
 
void slotForUpdate (const BeginOfRun *iT)
 
virtual ~Observer ()
 
- Public Member Functions inherited from Observer< const BeginOfEvent * >
 Observer ()
 
void slotForUpdate (const BeginOfEvent *iT)
 
virtual ~Observer ()
 
- Public Member Functions inherited from Observer< const BeginOfTrack * >
 Observer ()
 
void slotForUpdate (const BeginOfTrack *iT)
 
virtual ~Observer ()
 
- Public Member Functions inherited from Observer< const EndOfTrack * >
 Observer ()
 
void slotForUpdate (const EndOfTrack *iT)
 
virtual ~Observer ()
 
- Public Member Functions inherited from Observer< const EndOfEvent * >
 Observer ()
 
void slotForUpdate (const EndOfEvent *iT)
 
virtual ~Observer ()
 

Private Member Functions

double crystalDepth (G4LogicalVolume *, const G4ThreeVector &)
 
double crystalLength (G4LogicalVolume *)
 
double curve_LY (G4Step *)
 
void getBaseNumber (const G4Step *)
 
double getBirkL3 (G4Step *)
 
std::vector< double > getDDDArray (const std::string &, const DDsvalues_type &)
 
std::vector< std::string > getStringArray (const std::string &, const DDsvalues_type &)
 
void initMap (G4String, const DDCompactView &)
 

Private Attributes

EnergyResolutionVsLumi ageing
 
bool ageingWithSlopeLY
 
double birk1
 
double birk2
 
double birk3
 
double birkCut
 
double birkSlope
 
std::string crystalMat
 
std::string depth1Name
 
std::string depth2Name
 
bool ignoreDepCorr
 
std::vector< G4LogicalVolume * > noWeight
 
EcalNumberingSchemenumberingScheme
 
double slopeLY
 
bool storeRL
 
bool storeTrack
 
EcalBaseNumber theBaseNumber
 
bool useBirk
 
bool useBirkL3
 
std::vector< G4LogicalVolume * > useDepth1
 
std::vector< G4LogicalVolume * > useDepth2
 
bool useWeight
 
std::map< G4LogicalVolume
*, double > 
xtalLMap
 

Additional Inherited Members

- Public Types inherited from SensitiveDetector
enum  coordinates { WorldCoordinates, LocalCoordinates }
 
- Protected Member Functions inherited from CaloSD
G4bool checkHit ()
 
virtual void clearHits ()
 
CaloG4HitcreateNewHit ()
 
virtual bool filterHit (CaloG4Hit *, double)
 
double getAttenuation (G4Step *aStep, double birk1, double birk2, double birk3)
 
int getNumberOfHits ()
 
double getResponseWt (G4Track *)
 
virtual G4bool getStepInfo (G4Step *aStep)
 
G4bool hitExists ()
 
virtual void initRun ()
 
void resetForNewPrimary (const G4ThreeVector &, double)
 
G4ThreeVector setToGlobal (const G4ThreeVector &, const G4VTouchable *)
 
G4ThreeVector setToLocal (const G4ThreeVector &, const G4VTouchable *)
 
virtual void update (const BeginOfRun *)
 This routine will be called when the appropriate signal arrives. More...
 
virtual void update (const BeginOfEvent *)
 This routine will be called when the appropriate signal arrives. More...
 
virtual void update (const BeginOfTrack *trk)
 This routine will be called when the appropriate signal arrives. More...
 
virtual void update (const EndOfTrack *trk)
 This routine will be called when the appropriate signal arrives. More...
 
virtual void update (const ::EndOfEvent *)
 
void updateHit (CaloG4Hit *)
 
- Protected Member Functions inherited from Observer< const EndOfEvent * >
virtual void update (const EndOfEvent *)=0
 This routine will be called when the appropriate signal arrives. More...
 
- Protected Attributes inherited from CaloSD
int checkHits
 
double correctT
 
bool corrTOFBeam
 
CaloG4HitcurrentHit
 
CaloHitID currentID
 
float edepositEM
 
float edepositHAD
 
double eminHit
 
double eminHitD
 
G4int emPDG
 
double energyCut
 
G4ThreeVector entranceLocal
 
G4ThreeVector entrancePoint
 
G4int epPDG
 
bool forceSave
 
G4int gammaPDG
 
float incidentEnergy
 
double kmaxIon
 
double kmaxNeutron
 
double kmaxProton
 
const SimTrackManagerm_trackManager
 
G4ThreeVector posGlobal
 
G4StepPoint * preStepPoint
 
CaloHitID previousID
 
int primIDSaved
 
bool runInit
 
bool suppressHeavy
 
G4Track * theTrack
 
double tmaxHit
 
bool useMap
 

Detailed Description

Definition at line 26 of file ECalSD.h.

Constructor & Destructor Documentation

ECalSD::ECalSD ( G4String  name,
const DDCompactView cpv,
SensitiveDetectorCatalog clg,
edm::ParameterSet const &  p,
const SimTrackManager manager 
)

Definition at line 39 of file ECalSD.cc.

References DDFilteredView::addFilter(), ageing, ageingWithSlopeLY, birk1, birk2, birk3, birkCut, birkSlope, crystalMat, depth1Name, depth2Name, DDSpecificsFilter::equals, alcazmumu_cfi::filter, DDFilteredView::firstChild(), g, getDDDArray(), edm::ParameterSet::getParameter(), getStringArray(), edm::ParameterSet::getUntrackedParameter(), ignoreDepCorr, initMap(), CaloSD::kmaxIon, CaloSD::kmaxNeutron, CaloSD::kmaxProton, LogDebug, DDFilteredView::mergedSpecifics(), MeV, DDSpecificsFilter::setCriteria(), EnergyResolutionVsLumi::setLumies(), setNumberingScheme(), slopeLY, storeRL, storeTrack, AlCaHLTBitMon_QueryRunRegistry::string, CaloSD::suppressHeavy, useBirk, useBirkL3, and useWeight.

41  :
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")),
45  numberingScheme(0){
46 
47  // static SimpleConfigurable<bool> on1(false, "ECalSD:UseBirkLaw");
48  // static SimpleConfigurable<double> bk1(0.00463,"ECalSD:BirkC1");
49  // static SimpleConfigurable<double> bk2(-0.03, "ECalSD:BirkC2");
50  // static SimpleConfigurable<double> bk3(1.0, "ECalSD:BirkC3");
51  // Values from NIM A484 (2002) 239-244: as implemented in Geant3
52  // useBirk = on1.value();
53  // birk1 = bk1.value()*(g/(MeV*cm2));
54  // birk2 = bk2.value()*(g/(MeV*cm2))*(g/(MeV*cm2));
55 
56  edm::ParameterSet m_EC = p.getParameter<edm::ParameterSet>("ECalSD");
57  useBirk = m_EC.getParameter<bool>("UseBirkLaw");
58  useBirkL3 = m_EC.getParameter<bool>("BirkL3Parametrization");
59  birk1 = m_EC.getParameter<double>("BirkC1")*(g/(MeV*cm2));
60  birk2 = m_EC.getParameter<double>("BirkC2");
61  birk3 = m_EC.getParameter<double>("BirkC3");
62  birkSlope = m_EC.getParameter<double>("BirkSlope");
63  birkCut = m_EC.getParameter<double>("BirkCut");
64  slopeLY = m_EC.getParameter<double>("SlopeLightYield");
65  storeTrack = m_EC.getParameter<bool>("StoreSecondary");
66  ignoreDepCorr= m_EC.getParameter<bool>("IgnoreDepthCorr");
67  crystalMat = m_EC.getUntrackedParameter<std::string>("XtalMat","E_PbWO4");
68  bool isItTB = m_EC.getUntrackedParameter<bool>("TestBeam", false);
69  bool nullNS = m_EC.getUntrackedParameter<bool>("NullNumbering", false);
70  storeRL = m_EC.getUntrackedParameter<bool>("StoreRadLength", false);
71 
72  ageingWithSlopeLY = m_EC.getUntrackedParameter<bool>("AgeingWithSlopeLY", false);
73  if(ageingWithSlopeLY) ageing.setLumies(p.getParameter<edm::ParameterSet>("ECalSD").getParameter<double>("DelivLuminosity"),
74  p.getParameter<edm::ParameterSet>("ECalSD").getParameter<double>("InstLuminosity"));
75 
76  //Material list for HB/HE/HO sensitive detectors
77  std::string attribute = "ReadOutName";
79  DDValue ddv(attribute,name,0);
81  DDFilteredView fv(cpv);
82  fv.addFilter(filter);
83  fv.firstChild();
84  DDsvalues_type sv(fv.mergedSpecifics());
85  // Use of Weight
86  useWeight= true;
87  std::vector<double> tempD = getDDDArray("EnergyWeight",sv);
88  if (tempD.size() > 0) { if (tempD[0] < 0.1) useWeight = false; }
89  std::vector<std::string> tempS = getStringArray("Depth1Name",sv);
90  if (tempS.size() > 0) depth1Name = tempS[0];
91  else depth1Name = " ";
92  tempS = getStringArray("Depth2Name",sv);
93  if (tempS.size() > 0) depth2Name = tempS[0];
94  else depth2Name = " ";
95 
96  EcalNumberingScheme* scheme=0;
97  if (nullNS) scheme = 0;
98  else if (name == "EcalHitsEB") scheme = dynamic_cast<EcalNumberingScheme*>(new EcalBarrelNumberingScheme());
99  else if (name == "EcalHitsEE") scheme = dynamic_cast<EcalNumberingScheme*>(new EcalEndcapNumberingScheme());
100  else if (name == "EcalHitsES") {
101  if (isItTB) scheme = dynamic_cast<EcalNumberingScheme*>(new ESTBNumberingScheme());
102  else scheme = dynamic_cast<EcalNumberingScheme*>(new EcalPreshowerNumberingScheme());
103  useWeight = false;
104  } else {edm::LogWarning("EcalSim") << "ECalSD: ReadoutName not supported\n";}
105 
106  if (scheme) setNumberingScheme(scheme);
107 #ifdef DebugLog
108  LogDebug("EcalSim") << "Constructing a ECalSD with name " << GetName();
109 #endif
110  if (useWeight) {
111  edm::LogInfo("EcalSim") << "ECalSD:: Use of Birks law is set to "
112  << useBirk << " with three constants kB = "
113  << birk1 << ", C1 = " << birk2 << ", C2 = "
114  << birk3 <<"\n Use of L3 parametrization "
115  << useBirkL3 << " with slope " << birkSlope
116  << " and cut off " << birkCut << "\n"
117  << " Slope for Light yield is set to "
118  << slopeLY;
119  } else {
120  edm::LogInfo("EcalSim") << "ECalSD:: energy deposit is not corrected "
121  << " by Birk or light yield curve";
122  }
123 
124  edm::LogInfo("EcalSim") << "ECalSD:: Suppression Flag " << suppressHeavy
125  << " protons below " << kmaxProton << " MeV,"
126  << " neutrons below " << kmaxNeutron << " MeV and"
127  << " ions below " << kmaxIon << " MeV\n"
128  << " Depth1 Name = " << depth1Name
129  << " and Depth2 Name = " << depth2Name;
130 
131  if (useWeight) initMap(name,cpv);
132 
133 }
#define LogDebug(id)
bool ignoreDepCorr
Definition: ECalSD.h:55
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
bool useBirkL3
Definition: ECalSD.h:55
double kmaxNeutron
Definition: CaloSD.h:134
std::vector< std::string > getStringArray(const std::string &, const DDsvalues_type &)
Definition: ECalSD.cc:517
double birkSlope
Definition: ECalSD.h:56
double slopeLY
Definition: ECalSD.h:57
void setLumies(double x, double y)
double birk1
Definition: ECalSD.h:56
bool useWeight
Definition: ECalSD.h:54
double kmaxProton
Definition: CaloSD.h:134
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
Definition: Activities.doc:4
bool storeRL
Definition: ECalSD.h:54
CaloSD(G4String aSDname, const DDCompactView &cpv, SensitiveDetectorCatalog &clg, edm::ParameterSet const &p, const SimTrackManager *, int tSlice=1, bool ignoreTkID=false)
Definition: CaloSD.cc:24
const double MeV
std::string depth1Name
Definition: ECalSD.h:58
bool storeTrack
Definition: ECalSD.h:54
double kmaxIon
Definition: CaloSD.h:134
bool suppressHeavy
Definition: CaloSD.h:132
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::...
Definition: DDsvalues.h:19
void initMap(G4String, const DDCompactView &)
Definition: ECalSD.cc:277
std::string depth2Name
Definition: ECalSD.h:58
bool ageingWithSlopeLY
Definition: ECalSD.h:63
bool useBirk
Definition: ECalSD.h:55
std::vector< double > getDDDArray(const std::string &, const DDsvalues_type &)
Definition: ECalSD.cc:498
double birk2
Definition: ECalSD.h:56
void setCriteria(const DDValue &nameVal, comp_op, log_op l=AND, bool asString=true, bool merged=true)
Definition: DDFilter.cc:285
EcalNumberingScheme * numberingScheme
Definition: ECalSD.h:53
EnergyResolutionVsLumi ageing
Definition: ECalSD.h:62
void setNumberingScheme(EcalNumberingScheme *)
Definition: ECalSD.cc:267
double birkCut
Definition: ECalSD.h:56
std::string crystalMat
Definition: ECalSD.h:58
double birk3
Definition: ECalSD.h:56
The DDGenericFilter is a runtime-parametrized Filter looking on DDSpecifcs.
Definition: DDFilter.h:37
ECalSD::~ECalSD ( )
virtual

Definition at line 135 of file ECalSD.cc.

References numberingScheme.

135  {
136  if (numberingScheme) delete numberingScheme;
137 }
EcalNumberingScheme * numberingScheme
Definition: ECalSD.h:53

Member Function Documentation

double ECalSD::crystalDepth ( G4LogicalVolume *  lv,
const G4ThreeVector &  localPoint 
)
private

Definition at line 440 of file ECalSD.cc.

References funct::abs(), ignoreDepCorr, and xtalLMap.

Referenced by curve_LY(), and getRadiationLength().

441  {
442 
443  auto ite = xtalLMap.find(lv);
444  double depth(0);
445  if (ite != xtalLMap.end()) {
446  if (ignoreDepCorr) depth = (0.5*std::abs(ite->second)+localPoint.z());
447  else depth = std::abs(0.5*(ite->second)+localPoint.z());
448  }
449  return depth;
450 }
bool ignoreDepCorr
Definition: ECalSD.h:55
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::map< G4LogicalVolume *, double > xtalLMap
Definition: ECalSD.h:59
double ECalSD::crystalLength ( G4LogicalVolume *  lv)
private

Definition at line 434 of file ECalSD.cc.

References funct::abs(), and xtalLMap.

Referenced by curve_LY().

434  {
435 
436  auto ite = xtalLMap.find(lv);
437  double length = (ite == xtalLMap.end()) ? 230.0 : std::abs(ite->second);
438  return length;
439 }
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::map< G4LogicalVolume *, double > xtalLMap
Definition: ECalSD.h:59
double ECalSD::curve_LY ( G4Step *  aStep)
private

Definition at line 394 of file ECalSD.cc.

References ageing, ageingWithSlopeLY, EnergyResolutionVsLumi::calcLightCollectionEfficiencyWeighted(), crystalDepth(), crystalLength(), CaloSD::currentID, LogDebug, CaloSD::setToLocal(), slopeLY, CaloHitID::unitID(), and histoStyle::weight.

Referenced by getEnergyDeposit().

394  {
395 
396  G4StepPoint* stepPoint = aStep->GetPreStepPoint();
397  G4LogicalVolume* lv = stepPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
398 
399  double weight = 1.;
400  G4ThreeVector localPoint = setToLocal(stepPoint->GetPosition(),
401  stepPoint->GetTouchable());
402 
403  double crlength = crystalLength(lv);
404  double depth = crystalDepth(lv,localPoint);
405 
406  if(ageingWithSlopeLY){
407  //position along the crystal in mm from 0 to 230 (in EB)
408  if (depth >= -0.1 || depth <= crlength+0.1)
409  weight = ageing.calcLightCollectionEfficiencyWeighted(currentID.unitID(), depth/crlength);
410  }
411  else{
412  double dapd = crlength - depth;
413  if (dapd >= -0.1 || dapd <= crlength+0.1) {
414  if (dapd <= 100.)
415  weight = 1.0 + slopeLY - dapd * 0.01 * slopeLY;
416  } else {
417  edm::LogWarning("EcalSim") << "ECalSD: light coll curve : wrong distance "
418  << "to APD " << dapd << " crlength = "
419  << crlength <<" crystal name = " <<lv->GetName()
420  << " z of localPoint = " << localPoint.z()
421  << " take weight = " << weight;
422  }
423  }
424 #ifdef DebugLog
425  LogDebug("EcalSim") << "ECalSD, light coll curve : " << dapd
426  << " crlength = " << crlength
427  << " crystal name = " << lv->GetName()
428  << " z of localPoint = " << localPoint.z()
429  << " take weight = " << weight;
430 #endif
431  return weight;
432 }
#define LogDebug(id)
double slopeLY
Definition: ECalSD.h:57
double crystalLength(G4LogicalVolume *)
Definition: ECalSD.cc:434
bool ageingWithSlopeLY
Definition: ECalSD.h:63
double crystalDepth(G4LogicalVolume *, const G4ThreeVector &)
Definition: ECalSD.cc:440
CaloHitID currentID
Definition: CaloSD.h:117
double calcLightCollectionEfficiencyWeighted(DetId id, double z)
EnergyResolutionVsLumi ageing
Definition: ECalSD.h:62
int weight
Definition: histoStyle.py:50
uint32_t unitID() const
Definition: CaloHitID.h:22
G4ThreeVector setToLocal(const G4ThreeVector &, const G4VTouchable *)
Definition: CaloSD.cc:296
void ECalSD::getBaseNumber ( const G4Step *  aStep)
private

Definition at line 452 of file ECalSD.cc.

References EcalBaseNumber::addLevel(), EcalBaseNumber::getCapacity(), cuy::ii, LogDebug, EcalBaseNumber::reset(), EcalBaseNumber::setSize(), and theBaseNumber.

Referenced by setDetUnitId().

452  {
453 
455  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
456  int theSize = touch->GetHistoryDepth()+1;
457  if ( theBaseNumber.getCapacity() < theSize ) theBaseNumber.setSize(theSize);
458  //Get name and copy numbers
459  if ( theSize > 1 ) {
460  for (int ii = 0; ii < theSize ; ii++) {
461  theBaseNumber.addLevel(touch->GetVolume(ii)->GetName(),touch->GetReplicaNumber(ii));
462 #ifdef DebugLog
463  LogDebug("EcalSim") << "ECalSD::getBaseNumber(): Adding level " << ii
464  << ": " << touch->GetVolume(ii)->GetName() << "["
465  << touch->GetReplicaNumber(ii) << "]";
466 #endif
467  }
468  }
469 
470 }
#define LogDebug(id)
int ii
Definition: cuy.py:588
void addLevel(const std::string &name, const int &copyNumber)
EcalBaseNumber theBaseNumber
Definition: ECalSD.h:61
void setSize(const int &size)
double ECalSD::getBirkL3 ( G4Step *  aStep)
private

Definition at line 472 of file ECalSD.cc.

References birk1, birkCut, birkSlope, DeDxDiscriminatorTools::charge(), log, LogDebug, and histoStyle::weight.

Referenced by getEnergyDeposit().

472  {
473 
474  double weight = 1.;
475  double charge = aStep->GetPreStepPoint()->GetCharge();
476 
477  if (charge != 0. && aStep->GetStepLength() > 0) {
478  G4Material* mat = aStep->GetPreStepPoint()->GetMaterial();
479  double density = mat->GetDensity();
480  double dedx = aStep->GetTotalEnergyDeposit()/aStep->GetStepLength();
481  double rkb = birk1/density;
482  if (dedx > 0) {
483  weight = 1. - birkSlope*log(rkb*dedx);
484  if (weight < birkCut) weight = birkCut;
485  else if (weight > 1.) weight = 1.;
486  }
487 #ifdef DebugLog
488  LogDebug("EcalSim") << "ECalSD::getBirkL3 in " << mat->GetName()
489  << " Charge " << charge << " dE/dx " << dedx
490  << " Birk Const " << rkb << " Weight = " << weight
491  << " dE " << aStep->GetTotalEnergyDeposit();
492 #endif
493  }
494  return weight;
495 
496 }
#define LogDebug(id)
static std::vector< std::string > checklist log
double birkSlope
Definition: ECalSD.h:56
double birk1
Definition: ECalSD.h:56
double charge(const std::vector< uint8_t > &Ampls)
int weight
Definition: histoStyle.py:50
double birkCut
Definition: ECalSD.h:56
std::vector< double > ECalSD::getDDDArray ( const std::string &  str,
const DDsvalues_type sv 
)
private

Definition at line 498 of file ECalSD.cc.

References DDfetch(), DDValue::doubles(), LogDebug, and relativeConstraints::value.

Referenced by ECalSD().

499  {
500 
501 #ifdef DebugLog
502  LogDebug("EcalSim") << "ECalSD:getDDDArray called for " << str;
503 #endif
504  DDValue value(str);
505  if (DDfetch(&sv,value)) {
506 #ifdef DebugLog
507  LogDebug("EcalSim") << value;
508 #endif
509  const std::vector<double> & fvec = value.doubles();
510  return fvec;
511  } else {
512  std::vector<double> fvec;
513  return fvec;
514  }
515 }
#define LogDebug(id)
bool DDfetch(const DDsvalues_type *, DDValue &)
helper for retrieving DDValues from DDsvalues_type *.
Definition: DDsvalues.cc:102
uint16_t ECalSD::getDepth ( G4Step *  aStep)
virtual

Reimplemented from CaloSD.

Definition at line 227 of file ECalSD.cc.

References any(), getRadiationLength(), PCaloHit::kEcalDepthMask, PCaloHit::kEcalDepthOffset, PCaloHit::kEcalDepthRefz, storeRL, useDepth1, useDepth2, and xtalLMap.

227  {
228  G4LogicalVolume* lv = aStep->GetPreStepPoint()->GetTouchable()->GetVolume(0)->GetLogicalVolume();
229  uint16_t depth = any(useDepth1,lv) ? 1 : (any(useDepth2,lv) ? 2 : 0);
230  if (storeRL) {
231  auto ite = xtalLMap.find(lv);
232  uint16_t depth1 = (ite == xtalLMap.end()) ? 0 : (((ite->second) >= 0) ? 0 :
234  uint16_t depth2 = getRadiationLength(aStep);
235  depth |= (((depth2&PCaloHit::kEcalDepthMask) << PCaloHit::kEcalDepthOffset) | depth1);
236  }
237  return depth;
238 }
std::vector< G4LogicalVolume * > useDepth1
Definition: ECalSD.h:60
static const int kEcalDepthRefz
Definition: PCaloHit.h:70
bool any(const std::vector< T > &v, const T &what)
Definition: ECalSD.cc:34
bool storeRL
Definition: ECalSD.h:54
std::vector< G4LogicalVolume * > useDepth2
Definition: ECalSD.h:60
static const int kEcalDepthMask
Definition: PCaloHit.h:68
static const int kEcalDepthOffset
Definition: PCaloHit.h:69
virtual uint16_t getRadiationLength(G4Step *)
Definition: ECalSD.cc:240
std::map< G4LogicalVolume *, double > xtalLMap
Definition: ECalSD.h:59
double ECalSD::getEnergyDeposit ( G4Step *  aStep)
virtual

Reimplemented from CaloSD.

Definition at line 139 of file ECalSD.cc.

References any(), birk1, birk2, birk3, curve_LY(), CaloSD::getAttenuation(), getBirkL3(), CaloSD::getResponseWt(), TrackInformation::isPrimary(), ke, CaloSD::kmaxIon, CaloSD::kmaxNeutron, CaloSD::kmaxProton, LogDebug, MeV, noWeight, NULL, CaloSD::preStepPoint, CaloSD::suppressHeavy, CaloSD::theTrack, useBirk, useBirkL3, useWeight, and histoStyle::weight.

139  {
140 
141  if (aStep == NULL) {
142  return 0;
143  } else {
144  preStepPoint = aStep->GetPreStepPoint();
145  G4Track* theTrack = aStep->GetTrack();
146  double wt2 = theTrack->GetWeight();
147  G4String nameVolume = preStepPoint->GetPhysicalVolume()->GetName();
148 
149  // take into account light collection curve for crystals
150  double weight = 1.;
151  if (suppressHeavy) {
152  TrackInformation * trkInfo = (TrackInformation *)(theTrack->GetUserInformation());
153  if (trkInfo) {
154  int pdg = theTrack->GetDefinition()->GetPDGEncoding();
155  if (!(trkInfo->isPrimary())) { // Only secondary particles
156  double ke = theTrack->GetKineticEnergy()/MeV;
157  if (((pdg/1000000000 == 1 && ((pdg/10000)%100) > 0 &&
158  ((pdg/10)%100) > 0)) && (ke<kmaxIon)) weight = 0;
159  if ((pdg == 2212) && (ke < kmaxProton)) weight = 0;
160  if ((pdg == 2112) && (ke < kmaxNeutron)) weight = 0;
161 #ifdef DebugLog
162  if (weight == 0)
163  LogDebug("EcalSim") << "Ignore Track " << theTrack->GetTrackID()
164  << " Type " << theTrack->GetDefinition()->GetParticleName()
165  << " Kinetic Energy " << ke << " MeV";
166 #endif
167  }
168  }
169  }
170  G4LogicalVolume* lv = aStep->GetPreStepPoint()->GetTouchable()->GetVolume(0)->GetLogicalVolume();
171  if (useWeight && !any(noWeight,lv)) {
172  weight *= curve_LY(aStep);
173  if (useBirk) {
174  if (useBirkL3) weight *= getBirkL3(aStep);
175  else weight *= getAttenuation(aStep, birk1, birk2, birk3);
176  }
177  }
178  double wt1 = getResponseWt(theTrack);
179  double edep = aStep->GetTotalEnergyDeposit()*weight*wt1;
180  /*
181  if(wt2 != 1.0) {
182  std::cout << "ECalSD:: " << nameVolume
183  <<" LightWeight= " <<weight << " wt1= " <<wt1
184  << " wt2= " << wt2 << " "
185  << " Weighted Energy Deposit " << edep/MeV << " MeV"
186  << std::endl;
187  std::cout << theTrack->GetDefinition()->GetParticleName()
188  << " " << theTrack->GetKineticEnergy()
189  << " Id=" << theTrack->GetTrackID()
190  << " IdP=" << theTrack->GetParentID();
191  const G4VProcess* pr = theTrack->GetCreatorProcess();
192  if(pr) std::cout << " from " << pr->GetProcessName();
193  std::cout << std::endl;
194  }
195  */
196  if(wt2 > 0.0) { edep *= wt2; }
197 #ifdef DebugLog
198  LogDebug("EcalSim") << "ECalSD:: " << nameVolume
199  <<" Light Collection Efficiency " <<weight << ":" <<wt1
200  << " Weighted Energy Deposit " << edep/MeV << " MeV";
201 #endif
202  return edep;
203  }
204 }
#define LogDebug(id)
bool useBirkL3
Definition: ECalSD.h:55
double kmaxNeutron
Definition: CaloSD.h:134
bool any(const std::vector< T > &v, const T &what)
Definition: ECalSD.cc:34
double birk1
Definition: ECalSD.h:56
bool useWeight
Definition: ECalSD.h:54
#define NULL
Definition: scimark2.h:8
double kmaxProton
Definition: CaloSD.h:134
std::vector< G4LogicalVolume * > noWeight
Definition: ECalSD.h:60
const double MeV
double kmaxIon
Definition: CaloSD.h:134
bool suppressHeavy
Definition: CaloSD.h:132
double getAttenuation(G4Step *aStep, double birk1, double birk2, double birk3)
Definition: CaloSD.cc:465
G4Track * theTrack
Definition: CaloSD.h:118
double curve_LY(G4Step *)
Definition: ECalSD.cc:394
G4StepPoint * preStepPoint
Definition: CaloSD.h:120
bool useBirk
Definition: ECalSD.h:55
double getBirkL3(G4Step *)
Definition: ECalSD.cc:472
bool isPrimary() const
double getResponseWt(G4Track *)
Definition: CaloSD.cc:607
double birk2
Definition: ECalSD.h:56
int weight
Definition: histoStyle.py:50
double birk3
Definition: ECalSD.h:56
uint16_t ECalSD::getRadiationLength ( G4Step *  aStep)
virtual

Definition at line 240 of file ECalSD.cc.

References crystalDepth(), NULL, CaloSD::setToLocal(), and useWeight.

Referenced by getDepth().

240  {
241 
242  uint16_t thisX0 = 0;
243  if (aStep != NULL) {
244  G4StepPoint* hitPoint = aStep->GetPreStepPoint();
245  G4LogicalVolume* lv = hitPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
246 
247  if (useWeight) {
248  G4ThreeVector localPoint = setToLocal(hitPoint->GetPosition(),
249  hitPoint->GetTouchable());
250  double radl = hitPoint->GetMaterial()->GetRadlen();
251  double depth = crystalDepth(lv,localPoint);
252  thisX0 = (uint16_t)floor(depth/radl);
253  }
254  }
255  return thisX0;
256 }
bool useWeight
Definition: ECalSD.h:54
#define NULL
Definition: scimark2.h:8
double crystalDepth(G4LogicalVolume *, const G4ThreeVector &)
Definition: ECalSD.cc:440
G4ThreeVector setToLocal(const G4ThreeVector &, const G4VTouchable *)
Definition: CaloSD.cc:296
std::vector< std::string > ECalSD::getStringArray ( const std::string &  str,
const DDsvalues_type sv 
)
private

Definition at line 517 of file ECalSD.cc.

References DDfetch(), LogDebug, DDValue::strings(), and relativeConstraints::value.

Referenced by ECalSD().

518  {
519 
520 #ifdef DebugLog
521  LogDebug("EcalSim") << "ECalSD:getStringArray called for " << str;
522 #endif
523  DDValue value(str);
524  if (DDfetch(&sv,value)) {
525 #ifdef DebugLog
526  LogDebug("EcalSim") << value;
527 #endif
528  const std::vector<std::string> & fvec = value.strings();
529  return fvec;
530  } else {
531  std::vector<std::string> fvec;
532  return fvec;
533  }
534 }
#define LogDebug(id)
bool DDfetch(const DDsvalues_type *, DDValue &)
helper for retrieving DDValues from DDsvalues_type *.
Definition: DDsvalues.cc:102
int ECalSD::getTrackID ( G4Track *  aTrack)
virtual

Reimplemented from CaloSD.

Definition at line 206 of file ECalSD.cc.

References any(), archive::flag, CaloSD::forceSave, CaloSD::getTrackID(), CaloSD::preStepPoint, storeTrack, useDepth1, and useDepth2.

206  {
207 
208  int primaryID(0);
209  bool flag(false);
210  if (storeTrack) {
211  G4LogicalVolume* lv = preStepPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
212  if (any(useDepth1,lv)) {
213  flag = true;
214  } else if (any(useDepth2,lv)) {
215  flag = true;
216  }
217  }
218  if (flag) {
219  forceSave = true;
220  primaryID = aTrack->GetTrackID();
221  } else {
222  primaryID = CaloSD::getTrackID(aTrack);
223  }
224  return primaryID;
225 }
std::vector< G4LogicalVolume * > useDepth1
Definition: ECalSD.h:60
bool any(const std::vector< T > &v, const T &what)
Definition: ECalSD.cc:34
virtual int getTrackID(G4Track *)
Definition: CaloSD.cc:574
std::vector< G4LogicalVolume * > useDepth2
Definition: ECalSD.h:60
bool forceSave
Definition: CaloSD.h:137
bool storeTrack
Definition: ECalSD.h:54
G4StepPoint * preStepPoint
Definition: CaloSD.h:120
void ECalSD::initMap ( G4String  sd,
const DDCompactView cpv 
)
private

Definition at line 277 of file ECalSD.cc.

References DDFilteredView::addFilter(), any(), crystalMat, ddtrap, depth1Name, depth2Name, DDSpecificsFilter::equals, alcazmumu_cfi::filter, DDFilteredView::firstChild(), i, LogDebug, DDFilteredView::logicalPart(), DDLogicalPart::material(), DDName::name(), SensitiveDetector::name, DDBase< N, C >::name(), DDFilteredView::next(), noWeight, DDSolid::parameters(), DDSpecificsFilter::setCriteria(), DDSolid::shape(), DDLogicalPart::solid(), AlCaHLTBitMon_QueryRunRegistry::string, useDepth1, useDepth2, and xtalLMap.

Referenced by ECalSD().

277  {
278 
279  G4String attribute = "ReadOutName";
281  DDValue ddv(attribute,sd,0);
283  DDFilteredView fv(cpv);
284  fv.addFilter(filter);
285  fv.firstChild();
286 
287  std::vector<G4LogicalVolume*> lvused;
288  const G4LogicalVolumeStore * lvs = G4LogicalVolumeStore::GetInstance();
289  std::vector<G4LogicalVolume *>::const_iterator lvcite;
290  std::map<std::string, G4LogicalVolume *> nameMap;
291  for (auto lvi = lvs->begin(), lve = lvs->end(); lvi != lve; ++lvi)
292  nameMap.insert(std::make_pair((*lvi)->GetName(), *lvi));
293 
294  bool dodet=true;
295  while (dodet) {
296  const std::string &matname = fv.logicalPart().material().name().name();
297  const std::string &lvname = fv.logicalPart().name().name();
298  G4LogicalVolume* lv = nameMap[lvname];
299  int ibec = (lvname.find("EFRY") == std::string::npos) ? 0 : 1;
300  int iref = (lvname.find("refl") == std::string::npos) ? 0 : 1;
301  int type = (ibec+iref == 1) ? 1 : -1;
302  if (depth1Name != " ") {
303  if (strncmp(lvname.c_str(), depth1Name.c_str(), 4) == 0) {
304  if (!any(useDepth1, lv)) {
305  useDepth1.push_back(lv);
306 #ifdef DebugLog
307  LogDebug("EcalSim") << "ECalSD::initMap Logical Volume " << lvname
308  <<" in Depth 1 volume list";
309 #endif
310  }
311  G4LogicalVolume* lvr = nameMap[lvname + "_refl"];
312  if (lvr != 0 && !any(useDepth1, lvr)) {
313  useDepth1.push_back(lvr);
314 #ifdef DebugLog
315  LogDebug("EcalSim") << "ECalSD::initMap Logical Volume " << lvname << "_refl"
316  <<" in Depth 1 volume list";
317 #endif
318  }
319  }
320  }
321  if (depth2Name != " ") {
322  if (strncmp(lvname.c_str(), depth2Name.c_str(), 4) == 0) {
323  if (!any(useDepth2, lv)) {
324  useDepth2.push_back(lv);
325 #ifdef DebugLog
326  LogDebug("EcalSim") << "ECalSD::initMap Logical Volume " << lvname
327  <<" in Depth 2 volume list";
328 #endif
329  }
330  G4LogicalVolume* lvr = nameMap[lvname + "_refl"];
331  if (lvr != 0 && !any(useDepth2,lvr)) {
332  useDepth2.push_back(lvr);
333 #ifdef DebugLog
334  LogDebug("EcalSim") << "ECalSD::initMap Logical Volume " << lvname << "_refl"
335  <<" in Depth 2 volume list";
336 #endif
337  }
338  }
339  }
340  if (lv != 0) {
341  if (crystalMat.size() == matname.size() && !strcmp(crystalMat.c_str(), matname.c_str())) {
342  if (!any(lvused,lv)) {
343  lvused.push_back(lv);
344  const DDSolid & sol = fv.logicalPart().solid();
345  const std::vector<double> & paras = sol.parameters();
346 #ifdef DebugLog
347  LogDebug("EcalSim") << "ECalSD::initMap (for " << sd << "): Solid "
348  << lvname << " Shape " << sol.shape()
349  << " Parameter 0 = "<< paras[0]
350  << " Logical Volume " << lv;
351 #endif
352  if (sol.shape() == ddtrap) {
353  double dz = 2*paras[0];
354  xtalLMap.insert(std::pair<G4LogicalVolume*,double>(lv,dz*type));
355  lv = nameMap[lvname + "_refl"];
356  if (lv != 0)
357  xtalLMap.insert(std::pair<G4LogicalVolume*,double>(lv,-dz*type));
358  }
359  }
360  } else {
361  if (!any(noWeight,lv)) {
362  noWeight.push_back(lv);
363 #ifdef DebugLog
364  LogDebug("EcalSim") << "ECalSD::initMap Logical Volume " << lvname
365  << " Material " << matname <<" in noWeight list";
366 #endif
367  }
368  lv = nameMap[lvname];
369  if (lv != 0 && !any(noWeight,lv)) {
370  noWeight.push_back(lv);
371 #ifdef DebugLog
372  LogDebug("EcalSim") << "ECalSD::initMap Logical Volume " << lvname
373  << " Material " << matname <<" in noWeight list";
374 #endif
375  }
376  }
377  }
378  dodet = fv.next();
379  }
380 #ifdef DebugLog
381  LogDebug("EcalSim") << "ECalSD: Length Table for " << attribute << " = "
382  << sd << ":";
383  std::map<G4LogicalVolume*,double>::const_iterator ite = xtalLMap.begin();
384  int i=0;
385  for (; ite != xtalLMap.end(); ite++, i++) {
386  G4String name = "Unknown";
387  if (ite->first != 0) name = (ite->first)->GetName();
388  LogDebug("EcalSim") << " " << i << " " << ite->first << " " << name
389  << " L = " << ite->second;
390  }
391 #endif
392 }
#define LogDebug(id)
type
Definition: HCALResponse.h:21
int i
Definition: DBlmapReader.cc:9
const std::vector< double > & parameters(void) const
Give the parameters of the solid.
Definition: DDSolid.cc:150
std::vector< G4LogicalVolume * > useDepth1
Definition: ECalSD.h:60
bool any(const std::vector< T > &v, const T &what)
Definition: ECalSD.cc:34
std::vector< G4LogicalVolume * > noWeight
Definition: ECalSD.h:60
std::vector< G4LogicalVolume * > useDepth2
Definition: ECalSD.h:60
A DDSolid represents the shape of a part.
Definition: DDSolid.h:35
std::string depth1Name
Definition: ECalSD.h:58
DDSolidShape shape(void) const
The type of the solid.
Definition: DDSolid.cc:144
std::string depth2Name
Definition: ECalSD.h:58
double sd
void setCriteria(const DDValue &nameVal, comp_op, log_op l=AND, bool asString=true, bool merged=true)
Definition: DDFilter.cc:285
std::map< G4LogicalVolume *, double > xtalLMap
Definition: ECalSD.h:59
std::string crystalMat
Definition: ECalSD.h:58
The DDGenericFilter is a runtime-parametrized Filter looking on DDSpecifcs.
Definition: DDFilter.h:37
uint32_t ECalSD::setDetUnitId ( G4Step *  aStep)
virtual

Implements CaloSD.

Definition at line 258 of file ECalSD.cc.

References getBaseNumber(), EcalNumberingScheme::getUnitID(), numberingScheme, and theBaseNumber.

258  {
259  if (numberingScheme == 0) {
260  return EBDetId(1,1)();
261  } else {
262  getBaseNumber(aStep);
264  }
265 }
void getBaseNumber(const G4Step *)
Definition: ECalSD.cc:452
EcalBaseNumber theBaseNumber
Definition: ECalSD.h:61
EcalNumberingScheme * numberingScheme
Definition: ECalSD.h:53
virtual uint32_t getUnitID(const EcalBaseNumber &baseNumber) const =0
void ECalSD::setNumberingScheme ( EcalNumberingScheme scheme)

Definition at line 267 of file ECalSD.cc.

References numberingScheme.

Referenced by ECalSD(), and HcalTB04Analysis::update().

267  {
268  if (scheme != 0) {
269  edm::LogInfo("EcalSim") << "EcalSD: updates numbering scheme for "
270  << GetName() << "\n";
271  if (numberingScheme) delete numberingScheme;
272  numberingScheme = scheme;
273  }
274 }
EcalNumberingScheme * numberingScheme
Definition: ECalSD.h:53

Member Data Documentation

EnergyResolutionVsLumi ECalSD::ageing
private

Definition at line 62 of file ECalSD.h.

Referenced by curve_LY(), and ECalSD().

bool ECalSD::ageingWithSlopeLY
private

Definition at line 63 of file ECalSD.h.

Referenced by curve_LY(), and ECalSD().

double ECalSD::birk1
private

Definition at line 56 of file ECalSD.h.

Referenced by ECalSD(), getBirkL3(), and getEnergyDeposit().

double ECalSD::birk2
private

Definition at line 56 of file ECalSD.h.

Referenced by ECalSD(), and getEnergyDeposit().

double ECalSD::birk3
private

Definition at line 56 of file ECalSD.h.

Referenced by ECalSD(), and getEnergyDeposit().

double ECalSD::birkCut
private

Definition at line 56 of file ECalSD.h.

Referenced by ECalSD(), and getBirkL3().

double ECalSD::birkSlope
private

Definition at line 56 of file ECalSD.h.

Referenced by ECalSD(), and getBirkL3().

std::string ECalSD::crystalMat
private

Definition at line 58 of file ECalSD.h.

Referenced by ECalSD(), and initMap().

std::string ECalSD::depth1Name
private

Definition at line 58 of file ECalSD.h.

Referenced by ECalSD(), and initMap().

std::string ECalSD::depth2Name
private

Definition at line 58 of file ECalSD.h.

Referenced by ECalSD(), and initMap().

bool ECalSD::ignoreDepCorr
private

Definition at line 55 of file ECalSD.h.

Referenced by crystalDepth(), and ECalSD().

std::vector<G4LogicalVolume*> ECalSD::noWeight
private

Definition at line 60 of file ECalSD.h.

Referenced by getEnergyDeposit(), and initMap().

EcalNumberingScheme* ECalSD::numberingScheme
private

Definition at line 53 of file ECalSD.h.

Referenced by setDetUnitId(), setNumberingScheme(), and ~ECalSD().

double ECalSD::slopeLY
private

Definition at line 57 of file ECalSD.h.

Referenced by curve_LY(), and ECalSD().

bool ECalSD::storeRL
private

Definition at line 54 of file ECalSD.h.

Referenced by ECalSD(), and getDepth().

bool ECalSD::storeTrack
private

Definition at line 54 of file ECalSD.h.

Referenced by ECalSD(), and getTrackID().

EcalBaseNumber ECalSD::theBaseNumber
private

Definition at line 61 of file ECalSD.h.

Referenced by getBaseNumber(), and setDetUnitId().

bool ECalSD::useBirk
private

Definition at line 55 of file ECalSD.h.

Referenced by ECalSD(), and getEnergyDeposit().

bool ECalSD::useBirkL3
private

Definition at line 55 of file ECalSD.h.

Referenced by ECalSD(), and getEnergyDeposit().

std::vector<G4LogicalVolume*> ECalSD::useDepth1
private

Definition at line 60 of file ECalSD.h.

Referenced by getDepth(), getTrackID(), and initMap().

std::vector<G4LogicalVolume*> ECalSD::useDepth2
private

Definition at line 60 of file ECalSD.h.

Referenced by getDepth(), getTrackID(), and initMap().

bool ECalSD::useWeight
private

Definition at line 54 of file ECalSD.h.

Referenced by ECalSD(), getEnergyDeposit(), and getRadiationLength().

std::map<G4LogicalVolume*,double> ECalSD::xtalLMap
private

Definition at line 59 of file ECalSD.h.

Referenced by crystalDepth(), crystalLength(), getDepth(), and initMap().