test
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 "G4Material.hh"
20 #include "CLHEP/Units/GlobalSystemOfUnits.h"
21 
23  const SensitiveDetectorCatalog & clg,
24  edm::ParameterSet const & p,
25  const SimTrackManager* manager) :
26  CaloSD(name, cpv, clg, p, manager) {
27 
28  // Values from NIM 80 (1970) 239-244: as implemented in Geant3
29  edm::ParameterSet m_HC = p.getParameter<edm::ParameterSet>("HcalTB06BeamSD");
30  useBirk = m_HC.getParameter<bool>("UseBirkLaw");
31  birk1 = m_HC.getParameter<double>("BirkC1")*(g/(MeV*cm2));
32  birk2 = m_HC.getParameter<double>("BirkC2");
33  birk3 = m_HC.getParameter<double>("BirkC3");
34 
35  edm::LogInfo("HcalTB06BeamSD") << "HcalTB06BeamSD:: Use of Birks law is set to "
36  << useBirk << " with three constants kB = "
37  << birk1 << ", C1 = " <<birk2 << ", C2 = " <<birk3;
38 
39  std::string attribute, value;
40 
41  // Wire Chamber volume names
42  attribute = "Volume";
43  value = "WireChamber";
44  DDSpecificsFilter filter1;
45  DDValue ddv1(attribute,value,0);
46  filter1.setCriteria(ddv1,DDCompOp::equals);
47  DDFilteredView fv1(cpv);
48  fv1.addFilter(filter1);
49  wcNames = getNames(fv1);
50  edm::LogInfo("HcalTB06BeamSD")
51  << "HcalTB06BeamSD:: Names to be tested for "
52  << attribute << " = " << value << ": " << wcNames.size() << " paths";
53  for (unsigned int i=0; i<wcNames.size(); i++)
54  edm::LogInfo("HcalTB06BeamSD") << "HcalTB06BeamSD:: (" << i << ") "
55  << wcNames[i];
56 
57  //Material list for scintillator detector
58  attribute = "ReadOutName";
59  DDSpecificsFilter filter2;
60  DDValue ddv2(attribute,name,0);
61  filter2.setCriteria(ddv2,DDCompOp::equals);
62  DDFilteredView fv2(cpv);
63  fv2.addFilter(filter2);
64  bool dodet = fv2.firstChild();
65 
66  std::vector<G4String> matNames;
67  std::vector<int> nocc;
68  while (dodet) {
69  const DDLogicalPart & log = fv2.logicalPart();
70  matName = log.material().name().name();
71  bool notIn = true;
72  for (unsigned int i=0; i<matNames.size(); i++) {
73  if (matName == matNames[i]) {notIn = false; nocc[i]++;}
74  }
75  if (notIn) {
76  matNames.push_back(matName);
77  nocc.push_back(0);
78  }
79  dodet = fv2.next();
80  }
81  if (matNames.size() > 0) {
82  matName = matNames[0];
83  int occ = nocc[0];
84  for (unsigned int i = 0; i < matNames.size(); i++) {
85  if (nocc[i] > occ) {
86  occ = nocc[i];
87  matName = matNames[i];
88  }
89  }
90  } else {
91  matName = "Not Found";
92  }
93 
94  edm::LogInfo("HcalTB06BeamSD")
95  << "HcalTB06BeamSD: Material name for "
96  << attribute << " = " << name << ":" << matName;
97  matScin = G4Material::GetMaterial(matName);
98 }
99 
101 
102 double HcalTB06BeamSD::getEnergyDeposit(G4Step* aStep) {
103 
104  double destep = aStep->GetTotalEnergyDeposit();
105  double weight = 1;
106  if (useBirk && matScin == aStep->GetPreStepPoint()->GetMaterial()) {
107  weight *= getAttenuation(aStep, birk1, birk2, birk3);
108  }
109  LogDebug("HcalTB06BeamSD")
110  << "HcalTB06BeamSD: Detector "
111  << aStep->GetPreStepPoint()->GetTouchable()->GetVolume()->GetName()
112  << " weight " << weight;
113  return weight*destep;
114 }
115 
116 uint32_t HcalTB06BeamSD::setDetUnitId(G4Step * aStep) {
117 
118  G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
119  const G4VTouchable* touch = preStepPoint->GetTouchable();
120  G4String name = preStepPoint->GetPhysicalVolume()->GetName();
121 
122  int det = 1;
123  int lay = 0, x = 0, y = 0;
124  if (!isItWireChamber(name)) {
125  lay = (touch->GetReplicaNumber(0));
126  } else {
127  det = 2;
128  lay = (touch->GetReplicaNumber(1));
129  G4ThreeVector hitPoint = preStepPoint->GetPosition();
130  G4ThreeVector localPoint = setToLocal(hitPoint, touch);
131  x = (int)(localPoint.x()/(0.2*mm));
132  y = (int)(localPoint.y()/(0.2*mm));
133  }
134 
135  return packIndex (det, lay, x, y);
136 }
137 
138 uint32_t HcalTB06BeamSD::packIndex(int det, int lay, int x, int y) {
139 
140  int ix = 0, ixx = x;
141  if (x < 0) { ix = 1; ixx =-x;}
142  int iy = 0, iyy = y;
143  if (y < 0) { iy = 1; iyy =-y;}
144  uint32_t idx = (det&15)<<28; //bits 28-31
145  idx += (lay&127)<<21; //bits 21-27
146  idx += (iy&1)<<19; //bit 19
147  idx += (iyy&511)<<10; //bits 10-18
148  idx += (ix&1)<<9; //bit 9
149  idx += (ixx&511); //bits 0-8
150 
151  LogDebug("HcalTB06BeamSD") << "HcalTB06BeamSD: Detector " << det << " Layer "
152  << lay << " x " << x << " " << ix << " " << ixx
153  << " y " << y << " " << iy << " " << iyy << " ID "
154  << std::hex << idx << std::dec;
155  return idx;
156 }
157 
158 void HcalTB06BeamSD::unpackIndex(const uint32_t & idx, int& det, int& lay,
159  int& x, int& y) {
160 
161  det = (idx>>28)&15;
162  lay = (idx>>21)&127;
163  y = (idx>>10)&511; if (((idx>>19)&1) == 1) y = -y;
164  x = (idx)&511; if (((idx>>9)&1) == 1) x = -x;
165 
166 }
167 
168 std::vector<G4String> HcalTB06BeamSD::getNames(DDFilteredView& fv) {
169 
170  std::vector<G4String> tmp;
171  bool dodet = fv.firstChild();
172  while (dodet) {
173  const DDLogicalPart & log = fv.logicalPart();
174  G4String name = log.name().name();
175  bool ok = true;
176  for (unsigned int i=0; i<tmp.size(); i++)
177  if (name == tmp[i]) ok = false;
178  if (ok) tmp.push_back(name);
179  dodet = fv.next();
180  }
181  return tmp;
182 }
183 
185 
186  std::vector<G4String>::const_iterator it = wcNames.begin();
187  for (; it != wcNames.end(); it++)
188  if (name == *it) return true;
189  return false;
190 }
#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: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
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
const G4Material * matScin
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