CMS 3D CMS Logo

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 
17 
18 #include "G4Step.hh"
19 #include "G4Track.hh"
20 #include "G4Material.hh"
21 #include "CLHEP/Units/GlobalSystemOfUnits.h"
22 
23 HcalTB06BeamSD::HcalTB06BeamSD(const G4String& name, const DDCompactView & cpv,
24  const SensitiveDetectorCatalog & clg,
25  edm::ParameterSet const & p,
26  const SimTrackManager* manager) :
27  CaloSD(name, cpv, clg, p, manager) {
28 
29  // Values from NIM 80 (1970) 239-244: as implemented in Geant3
30  edm::ParameterSet m_HC = p.getParameter<edm::ParameterSet>("HcalTB06BeamSD");
31  useBirk = m_HC.getParameter<bool>("UseBirkLaw");
32  birk1 = m_HC.getParameter<double>("BirkC1")*(g/(MeV*cm2));
33  birk2 = m_HC.getParameter<double>("BirkC2");
34  birk3 = m_HC.getParameter<double>("BirkC3");
35 
36  edm::LogInfo("HcalTB06BeamSD") << "HcalTB06BeamSD:: Use of Birks law is set to "
37  << useBirk << " with three constants kB = "
38  << birk1 << ", C1 = " <<birk2 << ", C2 = " <<birk3;
39 
40  std::string attribute, value;
41 
42  // Wire Chamber volume names
43  attribute = "Volume";
44  value = "WireChamber";
45  DDSpecificsMatchesValueFilter filter1{DDValue(attribute,value,0)};
46  DDFilteredView fv1(cpv,filter1);
47  wcNames = getNames(fv1);
48  edm::LogInfo("HcalTB06BeamSD")
49  << "HcalTB06BeamSD:: Names to be tested for "
50  << attribute << " = " << value << ": " << wcNames.size() << " paths";
51  for (unsigned int i=0; i<wcNames.size(); i++)
52  edm::LogInfo("HcalTB06BeamSD") << "HcalTB06BeamSD:: (" << i << ") "
53  << wcNames[i];
54 
55  //Material list for scintillator detector
56  attribute = "ReadOutName";
57  DDSpecificsMatchesValueFilter filter2{DDValue(attribute,name,0)};
58  DDFilteredView fv2(cpv,filter2);
59  bool dodet = fv2.firstChild();
60 
61  std::vector<G4String> matNames;
62  std::vector<int> nocc;
63  while (dodet) {
64  const DDLogicalPart & log = fv2.logicalPart();
65  matName = log.material().name().name();
66  bool notIn = true;
67  for (unsigned int i=0; i<matNames.size(); i++) {
68  if (matName == matNames[i]) {notIn = false; nocc[i]++;}
69  }
70  if (notIn) {
71  matNames.push_back(matName);
72  nocc.push_back(0);
73  }
74  dodet = fv2.next();
75  }
76  if (!matNames.empty()) {
77  matName = matNames[0];
78  int occ = nocc[0];
79  for (unsigned int i = 0; i < matNames.size(); i++) {
80  if (nocc[i] > occ) {
81  occ = nocc[i];
82  matName = matNames[i];
83  }
84  }
85  } else {
86  matName = "Not Found";
87  }
88 
89  edm::LogInfo("HcalTB06BeamSD")
90  << "HcalTB06BeamSD: Material name for "
91  << attribute << " = " << name << ":" << matName;
92 }
93 
95 
96 double HcalTB06BeamSD::getEnergyDeposit(G4Step* aStep) {
97 
98  double destep = aStep->GetTotalEnergyDeposit();
99  double weight = 1;
100  if (useBirk && matName == aStep->GetPreStepPoint()->GetMaterial()->GetName()) {
101  weight *= getAttenuation(aStep, birk1, birk2, birk3);
102  }
103  LogDebug("HcalTB06BeamSD")
104  << "HcalTB06BeamSD: Detector "
105  << aStep->GetPreStepPoint()->GetTouchable()->GetVolume()->GetName()
106  << " weight " << weight;
107  return weight*destep;
108 }
109 
110 uint32_t HcalTB06BeamSD::setDetUnitId(G4Step * aStep) {
111 
112  G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
113  const G4VTouchable* touch = preStepPoint->GetTouchable();
114  G4String name = preStepPoint->GetPhysicalVolume()->GetName();
115 
116  int det = 1;
117  int lay = 0, x = 0, y = 0;
118  if (!isItWireChamber(name)) {
119  lay = (touch->GetReplicaNumber(0));
120  } else {
121  det = 2;
122  lay = (touch->GetReplicaNumber(1));
123  G4ThreeVector localPoint = setToLocal(preStepPoint->GetPosition(), touch);
124  x = (int)(localPoint.x()/(0.2*mm));
125  y = (int)(localPoint.y()/(0.2*mm));
126  }
127 
128  return HcalTestBeamNumbering::packIndex (det, lay, x, y);
129 }
130 
131 std::vector<G4String> HcalTB06BeamSD::getNames(DDFilteredView& fv) {
132 
133  std::vector<G4String> tmp;
134  bool dodet = fv.firstChild();
135  while (dodet) {
136  const DDLogicalPart & log = fv.logicalPart();
137  G4String name = log.name().name();
138  bool ok = true;
139  for (unsigned int i=0; i<tmp.size(); i++)
140  if (name == tmp[i]) ok = false;
141  if (ok) tmp.push_back(name);
142  dodet = fv.next();
143  }
144  return tmp;
145 }
146 
148 
149  std::vector<G4String>::const_iterator it = wcNames.begin();
150  for (; it != wcNames.end(); it++)
151  if (name == *it) return true;
152  return false;
153 }
#define LogDebug(id)
T getParameter(std::string const &) const
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
Definition: CaloSD.h:42
Definition: weight.py:1
type of data representation of DDCompactView
Definition: DDCompactView.h:90
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
static uint32_t packIndex(int det, int lay, int x, int y)
const double MeV
bool next()
set current node to the next node in the filtered tree
A DDLogicalPart aggregates information concerning material, solid and sensitveness ...
Definition: DDLogicalPart.h:92
HcalTB06BeamSD(const G4String &, const DDCompactView &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
double getAttenuation(G4Step *aStep, double birk1, double birk2, double birk3)
Definition: CaloSD.cc:465
virtual std::vector< std::string > getNames()
~HcalTB06BeamSD() override
G4StepPoint * preStepPoint
Definition: CaloSD.h:120
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)
double getEnergyDeposit(G4Step *) override
uint32_t setDetUnitId(G4Step *step) override
const std::string & name() const
Returns the name.
Definition: DDName.cc:90
const DDMaterial & material(void) const
Returns a reference object of the material this LogicalPart is made of.
G4ThreeVector setToLocal(const G4ThreeVector &, const G4VTouchable *)
Definition: CaloSD.cc:296