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