CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HcalTB06BeamSD.cc
Go to the documentation of this file.
1 // File: HcalTB06BeamSD.cc
3 // Description: Sensitive Detector class for beam counters in TB06 setup
5 
16 
17 #include "G4Step.hh"
18 #include "G4Track.hh"
19 #include "CLHEP/Units/GlobalSystemOfUnits.h"
20 
22  const SensitiveDetectorCatalog & clg,
23  edm::ParameterSet const & p,
24  const SimTrackManager* manager) :
25  CaloSD(name, cpv, clg, p, manager) {
26 
27  // Values from NIM 80 (1970) 239-244: as implemented in Geant3
28  edm::ParameterSet m_HC = p.getParameter<edm::ParameterSet>("HcalTB06BeamSD");
29  useBirk = m_HC.getParameter<bool>("UseBirkLaw");
30  birk1 = m_HC.getParameter<double>("BirkC1")*(g/(MeV*cm2));
31  birk2 = m_HC.getParameter<double>("BirkC2");
32  birk3 = m_HC.getParameter<double>("BirkC3");
33 
34  LogDebug("HcalTBSim") <<"***************************************************"
35  << "\n"
36  <<"* *"
37  << "\n"
38  << "* Constructing a HcalTB06BeamSD with name "
39  << name << "\n"
40  <<"* *"
41  << "\n"
42  <<"***************************************************";
43 
44  edm::LogInfo("HcalTBSim") << "HcalTB06BeamSD:: Use of Birks law is set to "
45  << useBirk << " with three constants kB = "
46  << birk1 << ", C1 = " <<birk2 << ", C2 = " <<birk3;
47 
48  std::string attribute, value;
49 
50  // Wire Chamber volume names
51  attribute = "Volume";
52  value = "WireChamber";
53  DDSpecificsFilter filter1;
54  DDValue ddv1(attribute,value,0);
55  filter1.setCriteria(ddv1,DDCompOp::equals);
56  DDFilteredView fv1(cpv);
57  fv1.addFilter(filter1);
58  wcNames = getNames(fv1);
59  edm::LogInfo("HcalTBSim") << "HcalTB06BeamSD:: Names to be tested for "
60  << attribute << " = " << value << ": "
61  << wcNames.size() << " paths";
62  for (unsigned int i=0; i<wcNames.size(); i++)
63  edm::LogInfo("HcalTBSim") << "HcalTB06BeamSD:: (" << i << ") "
64  << wcNames[i];
65 
66  //Material list for scintillator detector
67  attribute = "ReadOutName";
68  DDSpecificsFilter filter2;
69  DDValue ddv2(attribute,name,0);
70  filter2.setCriteria(ddv2,DDCompOp::equals);
71  DDFilteredView fv2(cpv);
72  fv2.addFilter(filter2);
73  bool dodet = fv2.firstChild();
74 
75  std::vector<G4String> matNames;
76  std::vector<int> nocc;
77  while (dodet) {
78  const DDLogicalPart & log = fv2.logicalPart();
79  matName = log.material().name().name();
80  bool notIn = true;
81  for (unsigned int i=0; i<matNames.size(); i++)
82  if (matName == matNames[i]) {notIn = false; nocc[i]++;}
83  if (notIn) {
84  matNames.push_back(matName);
85  nocc.push_back(0);
86  }
87  dodet = fv2.next();
88  }
89  if (matNames.size() > 0) {
90  matName = matNames[0];
91  int occ = nocc[0];
92  for (unsigned int i = 0; i < matNames.size(); i++) {
93  if (nocc[i] > occ) {
94  occ = nocc[i];
95  matName = matNames[i];
96  }
97  }
98  } else {
99  matName = "Not Found";
100  }
101 
102  edm::LogInfo("HcalTBSim") << "HcalTB06BeamSD: Material name for "
103  << attribute << " = " << name << ":" << matName;
104 }
105 
107 
108 double HcalTB06BeamSD::getEnergyDeposit(G4Step* aStep) {
109 
110  double destep = aStep->GetTotalEnergyDeposit();
111  double weight = 1;
112  if (useBirk) {
113  G4Material* mat = aStep->GetPreStepPoint()->GetMaterial();
114  if (mat->GetName() == matName)
115  weight *= getAttenuation(aStep, birk1, birk2, birk3);
116  }
117  LogDebug("HcalTBSim") << "HcalTB06BeamSD: Detector "
118  << aStep->GetPreStepPoint()->GetTouchable()->GetVolume()->GetName()
119  << " weight " << weight;
120  return weight*destep;
121 }
122 
123 uint32_t HcalTB06BeamSD::setDetUnitId(G4Step * aStep) {
124 
125  G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
126  const G4VTouchable* touch = preStepPoint->GetTouchable();
127  G4String name = preStepPoint->GetPhysicalVolume()->GetName();
128 
129  int det = 1;
130  int lay = 0, x = 0, y = 0;
131  if (!isItWireChamber(name)) {
132  lay = (touch->GetReplicaNumber(0));
133  } else {
134  det = 2;
135  lay = (touch->GetReplicaNumber(1));
136  G4ThreeVector hitPoint = preStepPoint->GetPosition();
137  G4ThreeVector localPoint = setToLocal(hitPoint, touch);
138  x = (int)(localPoint.x()/(0.2*mm));
139  y = (int)(localPoint.y()/(0.2*mm));
140  }
141 
142  return packIndex (det, lay, x, y);
143 }
144 
145 uint32_t HcalTB06BeamSD::packIndex(int det, int lay, int x, int y) {
146 
147  int ix = 0, ixx = x;
148  if (x < 0) { ix = 1; ixx =-x;}
149  int iy = 0, iyy = y;
150  if (y < 0) { iy = 1; iyy =-y;}
151  uint32_t idx = (det&15)<<28; //bits 28-31
152  idx += (lay&127)<<21; //bits 21-27
153  idx += (iy&1)<<19; //bit 19
154  idx += (iyy&511)<<10; //bits 10-18
155  idx += (ix&1)<<9; //bit 9
156  idx += (ixx&511); //bits 0-8
157 
158  LogDebug("HcalTBSim") << "HcalTB06BeamSD: Detector " << det << " Layer "
159  << lay << " x " << x << " " << ix << " " << ixx
160  << " y " << y << " " << iy << " " << iyy << " ID "
161  << std::hex << idx << std::dec;
162  return idx;
163 }
164 
165 void HcalTB06BeamSD::unpackIndex(const uint32_t & idx, int& det, int& lay,
166  int& x, int& y) {
167 
168  det = (idx>>28)&15;
169  lay = (idx>>21)&127;
170  y = (idx>>10)&511; if (((idx>>19)&1) == 1) y = -y;
171  x = (idx)&511; if (((idx>>9)&1) == 1) x = -x;
172 
173 }
174 
175 std::vector<G4String> HcalTB06BeamSD::getNames(DDFilteredView& fv) {
176 
177  std::vector<G4String> tmp;
178  bool dodet = fv.firstChild();
179  while (dodet) {
180  const DDLogicalPart & log = fv.logicalPart();
181  G4String name = log.name().name();
182  bool ok = true;
183  for (unsigned int i=0; i<tmp.size(); i++)
184  if (name == tmp[i]) ok = false;
185  if (ok) tmp.push_back(name);
186  dodet = fv.next();
187  }
188  return tmp;
189 }
190 
192 
193  std::vector<G4String>::const_iterator it = wcNames.begin();
194  for (; it != wcNames.end(); it++)
195  if (name == *it) return true;
196  return false;
197 }
#define LogDebug(id)
virtual double getEnergyDeposit(G4Step *)
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
const DDLogicalPart & logicalPart() const
The logical-part of the current node in the filtered-view.
std::vector< G4String > wcNames
const N & name() const
Definition: DDBase.h:78
HcalTB06BeamSD(G4String, const DDCompactView &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
Definition: CaloSD.h:42
void addFilter(const DDFilter &, DDLogOp op=DDLogOp::AND)
int iyy[18][41][3]
int ixx[18][41][3]
type of data representation of DDCompactView
Definition: DDCompactView.h:76
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
virtual ~HcalTB06BeamSD()
const double MeV
T x() const
Cartesian x coordinate.
bool next()
set current node to the next node in the filtered tree
virtual uint32_t setDetUnitId(G4Step *step)
A DDLogicalPart aggregates information concerning material, solid and sensitveness ...
Definition: DDLogicalPart.h:88
double getAttenuation(G4Step *aStep, double birk1, double birk2, double birk3)
Definition: CaloSD.cc:465
virtual std::vector< std::string > getNames()
static uint32_t packIndex(int det, int lay, int x, int y)
G4StepPoint * preStepPoint
Definition: CaloSD.h:120
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
G4String matName
bool firstChild()
set the current node to the first child ...
bool isItWireChamber(G4String)
int weight
Definition: histoStyle.py:50
static void unpackIndex(const uint32_t &idx, int &det, int &lay, int &x, int &y)
const std::string & name() const
Returns the name.
Definition: DDName.cc:87
const DDMaterial & material(void) const
Returns a reference object of the material this LogicalPart is made of.
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