CMS 3D CMS Logo

HcalTB02SD.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: HcalTestBeam
4 // Class : HcalTB02SD
5 //
6 // Implementation:
7 // Sensitive Detector class for Hcal Test Beam 2002 detectors
8 //
9 // Original Author:
10 // Created: Sun 21 10:14:34 CEST 2006
11 //
12 
13 // system include files
14 
15 // user include files
25 
26 #include "G4Step.hh"
27 #include "G4Track.hh"
28 #include "G4VProcess.hh"
29 
30 #include "G4SystemOfUnits.hh"
31 
32 //
33 // constructors and destructor
34 //
35 
37  const SensitiveDetectorCatalog & clg,
38  edm::ParameterSet const & p,
39  const SimTrackManager* manager) :
40  CaloSD(name, cpv, clg, p, manager), numberingScheme(nullptr) {
41 
42  edm::ParameterSet m_SD = p.getParameter<edm::ParameterSet>("HcalTB02SD");
43  useBirk= m_SD.getUntrackedParameter<bool>("UseBirkLaw",false);
44  birk1 = m_SD.getUntrackedParameter<double>("BirkC1",0.013)*(g/(MeV*cm2));
45  birk2 = m_SD.getUntrackedParameter<double>("BirkC2",0.0568);
46  birk3 = m_SD.getUntrackedParameter<double>("BirkC3",1.75);
47  useWeight= true;
48 
49  HcalTB02NumberingScheme* scheme=nullptr;
50  if (name == "EcalHitsEB") {
51  scheme = dynamic_cast<HcalTB02NumberingScheme*>(new HcalTB02XtalNumberingScheme());
52  useBirk = false;
53  } else if (name == "HcalHits") {
54  scheme = dynamic_cast<HcalTB02NumberingScheme*>(new HcalTB02HcalNumberingScheme());
55  useWeight= false;
56  } else {edm::LogWarning("HcalTBSim") << "HcalTB02SD: ReadoutName " << name
57  << " not supported\n";}
58 
59  if (scheme) setNumberingScheme(scheme);
60  LogDebug("HcalTBSim")
61  << "***************************************************"
62  << "\n"
63  << "* *"
64  << "\n"
65  << "* Constructing a HcalTB02SD with name " << GetName()
66  << "\n"
67  << "* *"
68  << "\n"
69  << "***************************************************" ;
70  edm::LogInfo("HcalTBSim") << "HcalTB02SD:: Use of Birks law is set to "
71  << useBirk << " with three constants kB = "
72  << birk1 << ", C1 = " << birk2 << ", C2 = "
73  << birk3;
74 
75  if (useWeight) initMap(name,cpv);
76 
77 }
78 
80  if (numberingScheme) delete numberingScheme;
81 }
82 
83 //
84 // member functions
85 //
86 
87 double HcalTB02SD::getEnergyDeposit(const G4Step * aStep) {
88 
89  if (aStep == nullptr) {
90  return 0;
91  } else {
92  auto const preStepPoint = aStep->GetPreStepPoint();
93  auto const & nameVolume = preStepPoint->GetPhysicalVolume()->GetName();
94 
95  // take into account light collection curve for crystals
96  double weight = 1.;
97  if (useWeight) weight *= curve_LY(nameVolume, preStepPoint);
98  if (useBirk) weight *= getAttenuation(aStep, birk1, birk2, birk3);
99  double edep = aStep->GetTotalEnergyDeposit() * weight;
100  LogDebug("HcalTBSim") << "HcalTB02SD:: " << nameVolume
101  <<" Light Collection Efficiency " << weight
102  << " Weighted Energy Deposit " << edep/MeV << " MeV";
103  return edep;
104  }
105 }
106 
107 uint32_t HcalTB02SD::setDetUnitId(const G4Step * aStep) {
108  return (numberingScheme == nullptr ? 0 : (uint32_t)(numberingScheme->getUnitID(aStep)));
109 }
110 
112  if (scheme != nullptr) {
113  edm::LogInfo("HcalTBSim") << "HcalTB02SD: updates numbering scheme for "
114  << GetName();
115  if (numberingScheme) delete numberingScheme;
116  numberingScheme = scheme;
117  }
118 }
119 
120 void HcalTB02SD::initMap(const std::string& sd, const DDCompactView & cpv) {
121 
122  G4String attribute = "ReadOutName";
124  DDFilteredView fv(cpv,filter);
125  fv.firstChild();
126 
127  bool dodet=true;
128  while (dodet) {
129  const DDSolid & sol = fv.logicalPart().solid();
130  const std::vector<double> & paras = sol.parameters();
131  G4String name = sol.name().name();
132  LogDebug("HcalTBSim") << "HcalTB02SD::initMap (for " << sd << "): Solid "
133  << name << " Shape " << sol.shape()
134  << " Parameter 0 = " << paras[0];
135  if (sol.shape() == DDSolidShape::ddtrap) {
136  double dz = 2*paras[0];
137  lengthMap.insert(std::pair<G4String,double>(name,dz));
138  }
139  dodet = fv.next();
140  }
141  LogDebug("HcalTBSim") << "HcalTB02SD: Length Table for " << attribute
142  << " = " << sd << ":";
143  std::map<G4String,double>::const_iterator it = lengthMap.begin();
144  int i=0;
145  for (; it != lengthMap.end(); it++, i++) {
146  LogDebug("HcalTBSim") << " " << i << " " << it->first << " L = "
147  << it->second;
148  }
149 }
150 
151 double HcalTB02SD::curve_LY(const G4String& nameVolume, const G4StepPoint* stepPoint) {
152 
153  double weight = 1.;
154  G4ThreeVector localPoint = setToLocal(stepPoint->GetPosition(),
155  stepPoint->GetTouchable());
156  double crlength = crystalLength(nameVolume);
157  double dapd = 0.5 * crlength - localPoint.z();
158  if (dapd >= -0.1 || dapd <= crlength+0.1) {
159  if (dapd <= 100.)
160  weight = 1.05 - dapd * 0.0005;
161  } else {
162  edm::LogWarning("HcalTBSim") << "HcalTB02SD: light coll curve : wrong "
163  << "distance to APD " << dapd <<" crlength = "
164  << crlength << " crystal name = " <<nameVolume
165  << " z of localPoint = " << localPoint.z()
166  << " take weight = " << weight;
167  }
168  LogDebug("HcalTBSim") << "HcalTB02SD, light coll curve : " << dapd
169  << " crlength = " << crlength
170  << " crystal name = " << nameVolume
171  << " z of localPoint = " << localPoint.z()
172  << " take weight = " << weight;
173  return weight;
174 }
175 
176 double HcalTB02SD::crystalLength(const G4String& name) {
177 
178  double length = 230.;
179  std::map<G4String,double>::const_iterator it = lengthMap.find(name);
180  if (it != lengthMap.end()) length = it->second;
181  return length;
182 }
#define LogDebug(id)
std::map< G4String, double > lengthMap
Definition: HcalTB02SD.h:54
T getParameter(std::string const &) const
double birk2
Definition: HcalTB02SD.h:53
T getUntrackedParameter(std::string const &, T const &) const
const std::vector< double > & parameters(void) const
Give the parameters of the solid.
Definition: DDSolid.cc:144
HcalTB02NumberingScheme * numberingScheme
Definition: HcalTB02SD.h:50
const DDLogicalPart & logicalPart() const
The logical-part of the current node in the filtered-view.
const N & name() const
Definition: DDBase.h:74
double crystalLength(const G4String &)
Definition: HcalTB02SD.cc:176
Definition: CaloSD.h:37
double getEnergyDeposit(const G4Step *) override
Definition: HcalTB02SD.cc:87
Definition: weight.py:1
uint32_t setDetUnitId(const G4Step *step) override
Definition: HcalTB02SD.cc:107
bool useWeight
Definition: HcalTB02SD.h:51
const DDSolid & solid(void) const
Returns a reference object of the solid being the shape of this LogicalPart.
bool useBirk
Definition: HcalTB02SD.h:52
#define nullptr
Compact representation of the geometrical detector hierarchy.
Definition: DDCompactView.h:80
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
A DDSolid represents the shape of a part.
Definition: DDSolid.h:39
static TrackerG4SimHitNumberingScheme & numberingScheme(const DDCompactView &cpv, const GeometricDet &det)
~HcalTB02SD() override
Definition: HcalTB02SD.cc:79
const double MeV
double curve_LY(const G4String &, const G4StepPoint *)
Definition: HcalTB02SD.cc:151
double birk3
Definition: HcalTB02SD.h:53
bool next()
set current node to the next node in the filtered tree
DDSolidShape shape(void) const
The type of the solid.
Definition: DDSolid.cc:138
G4ThreeVector setToLocal(const G4ThreeVector &, const G4VTouchable *) const
Definition: CaloSD.cc:307
void setNumberingScheme(HcalTB02NumberingScheme *scheme)
Definition: HcalTB02SD.cc:111
void initMap(const std::string &, const DDCompactView &)
Definition: HcalTB02SD.cc:120
virtual int getUnitID(const G4Step *aStep) const =0
double sd
double getAttenuation(const G4Step *aStep, double birk1, double birk2, double birk3) const
Definition: CaloSD.cc:462
HcalTB02SD(const std::string &, const DDCompactView &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
Definition: HcalTB02SD.cc:36
bool firstChild()
set the current node to the first child ...
const std::string & name() const
Returns the name.
Definition: DDName.cc:53
double birk1
Definition: HcalTB02SD.h:53