CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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(0) {
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=0;
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(G4Step * aStep) {
88 
89  if (aStep == NULL) {
90  return 0;
91  } else {
92  preStepPoint = aStep->GetPreStepPoint();
93  G4String 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(G4Step * aStep) {
108  return (numberingScheme == 0 ? 0 : (uint32_t)(numberingScheme->getUnitID(aStep)));
109 }
110 
112  if (scheme != 0) {
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(G4String sd, const DDCompactView & cpv) {
121 
122  G4String attribute = "ReadOutName";
124  DDValue ddv(attribute,sd,0);
125  filter.setCriteria(ddv,DDCompOp::equals);
126  DDFilteredView fv(cpv);
127  fv.addFilter(filter);
128  fv.firstChild();
129 
130  bool dodet=true;
131  while (dodet) {
132  const DDSolid & sol = fv.logicalPart().solid();
133  const std::vector<double> & paras = sol.parameters();
134  G4String name = sol.name().name();
135  LogDebug("HcalTBSim") << "HcalTB02SD::initMap (for " << sd << "): Solid "
136  << name << " Shape " << sol.shape()
137  << " Parameter 0 = " << paras[0];
138  if (sol.shape() == ddtrap) {
139  double dz = 2*paras[0];
140  lengthMap.insert(std::pair<G4String,double>(name,dz));
141  }
142  dodet = fv.next();
143  }
144  LogDebug("HcalTBSim") << "HcalTB02SD: Length Table for " << attribute
145  << " = " << sd << ":";
146  std::map<G4String,double>::const_iterator it = lengthMap.begin();
147  int i=0;
148  for (; it != lengthMap.end(); it++, i++) {
149  LogDebug("HcalTBSim") << " " << i << " " << it->first << " L = "
150  << it->second;
151  }
152 }
153 
154 double HcalTB02SD::curve_LY(G4String& nameVolume, G4StepPoint* stepPoint) {
155 
156  double weight = 1.;
157  G4ThreeVector localPoint = setToLocal(stepPoint->GetPosition(),
158  stepPoint->GetTouchable());
159  double crlength = crystalLength(nameVolume);
160  double dapd = 0.5 * crlength - localPoint.z();
161  if (dapd >= -0.1 || dapd <= crlength+0.1) {
162  if (dapd <= 100.)
163  weight = 1.05 - dapd * 0.0005;
164  } else {
165  edm::LogWarning("HcalTBSim") << "HcalTB02SD: light coll curve : wrong "
166  << "distance to APD " << dapd <<" crlength = "
167  << crlength << " crystal name = " <<nameVolume
168  << " z of localPoint = " << localPoint.z()
169  << " take weight = " << weight;
170  }
171  LogDebug("HcalTBSim") << "HcalTB02SD, light coll curve : " << dapd
172  << " crlength = " << crlength
173  << " crystal name = " << nameVolume
174  << " z of localPoint = " << localPoint.z()
175  << " take weight = " << weight;
176  return weight;
177 }
178 
179 double HcalTB02SD::crystalLength(G4String name) {
180 
181  double length = 230.;
182  std::map<G4String,double>::const_iterator it = lengthMap.find(name);
183  if (it != lengthMap.end()) length = it->second;
184  return length;
185 }
#define LogDebug(id)
std::map< G4String, double > lengthMap
Definition: HcalTB02SD.h:51
T getParameter(std::string const &) const
double birk2
Definition: HcalTB02SD.h:50
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
const std::vector< double > & parameters(void) const
Give the parameters of the solid.
Definition: DDSolid.cc:150
HcalTB02NumberingScheme * numberingScheme
Definition: HcalTB02SD.h:47
const DDLogicalPart & logicalPart() const
The logical-part of the current node in the filtered-view.
const N & name() const
Definition: DDBase.h:78
virtual ~HcalTB02SD()
Definition: HcalTB02SD.cc:79
Definition: CaloSD.h:42
void addFilter(const DDFilter &, DDLogOp op=DDLogOp::AND)
#define NULL
Definition: scimark2.h:8
bool useWeight
Definition: HcalTB02SD.h:48
const DDSolid & solid(void) const
Returns a reference object of the solid being the shape of this LogicalPart.
bool useBirk
Definition: HcalTB02SD.h:49
type of data representation of DDCompactView
Definition: DDCompactView.h:77
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:35
static TrackerG4SimHitNumberingScheme & numberingScheme(const DDCompactView &cpv, const GeometricDet &det)
const double MeV
double birk3
Definition: HcalTB02SD.h:50
bool next()
set current node to the next node in the filtered tree
HcalTB02SD(G4String, const DDCompactView &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
Definition: HcalTB02SD.cc:36
DDSolidShape shape(void) const
The type of the solid.
Definition: DDSolid.cc:144
double crystalLength(G4String)
Definition: HcalTB02SD.cc:179
void initMap(G4String, const DDCompactView &)
Definition: HcalTB02SD.cc:120
double getAttenuation(G4Step *aStep, double birk1, double birk2, double birk3)
Definition: CaloSD.cc:465
virtual uint32_t setDetUnitId(G4Step *step)
Definition: HcalTB02SD.cc:107
G4StepPoint * preStepPoint
Definition: CaloSD.h:120
void setNumberingScheme(HcalTB02NumberingScheme *scheme)
Definition: HcalTB02SD.cc:111
double sd
bool firstChild()
set the current node to the first child ...
virtual int getUnitID(const G4Step *aStep) const =0
double curve_LY(G4String &, G4StepPoint *)
Definition: HcalTB02SD.cc:154
const std::string & name() const
Returns the name.
Definition: DDName.cc:87
virtual double getEnergyDeposit(G4Step *)
Definition: HcalTB02SD.cc:87
double birk1
Definition: HcalTB02SD.h:50
void setCriteria(const DDValue &nameVal, DDCompOp, DDLogOp l=DDLogOp::AND, bool asString=true, bool merged=true)
Definition: DDFilter.cc:245
The DDGenericFilter is a runtime-parametrized Filter looking on DDSpecifcs.
Definition: DDFilter.h:32
G4ThreeVector setToLocal(const G4ThreeVector &, const G4VTouchable *)
Definition: CaloSD.cc:296