00001
00002
00003
00005
00006 #include "SimG4CMS/HcalTestBeam/interface/HcalTB06BeamSD.h"
00007 #include "SimG4Core/Notification/interface/TrackInformation.h"
00008 #include "DetectorDescription/Core/interface/DDFilter.h"
00009 #include "DetectorDescription/Core/interface/DDFilteredView.h"
00010 #include "DetectorDescription/Core/interface/DDLogicalPart.h"
00011 #include "DetectorDescription/Core/interface/DDMaterial.h"
00012 #include "DetectorDescription/Core/interface/DDSplit.h"
00013 #include "DetectorDescription/Core/interface/DDValue.h"
00014 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00015 #include "FWCore/Utilities/interface/Exception.h"
00016
00017 #include "G4Step.hh"
00018 #include "G4Track.hh"
00019 #include "CLHEP/Units/SystemOfUnits.h"
00020
00021 HcalTB06BeamSD::HcalTB06BeamSD(G4String name, const DDCompactView & cpv,
00022 SensitiveDetectorCatalog & clg,
00023 edm::ParameterSet const & p,
00024 const SimTrackManager* manager) :
00025 CaloSD(name, cpv, clg, p, manager) {
00026
00027
00028 edm::ParameterSet m_HC = p.getParameter<edm::ParameterSet>("HcalTB06BeamSD");
00029 useBirk = m_HC.getParameter<bool>("UseBirkLaw");
00030 birk1 = m_HC.getParameter<double>("BirkC1")*(g/(MeV*cm2));
00031 birk2 = m_HC.getParameter<double>("BirkC2");
00032 birk3 = m_HC.getParameter<double>("BirkC3");
00033
00034 LogDebug("HcalTBSim") <<"***************************************************"
00035 << "\n"
00036 <<"* *"
00037 << "\n"
00038 << "* Constructing a HcalTB06BeamSD with name "
00039 << name << "\n"
00040 <<"* *"
00041 << "\n"
00042 <<"***************************************************";
00043
00044 edm::LogInfo("HcalTBSim") << "HcalTB06BeamSD:: Use of Birks law is set to "
00045 << useBirk << " with three constants kB = "
00046 << birk1 << ", C1 = " <<birk2 << ", C2 = " <<birk3;
00047
00048 std::string attribute, value;
00049
00050
00051 attribute = "Volume";
00052 value = "WireChamber";
00053 DDSpecificsFilter filter1;
00054 DDValue ddv1(attribute,value,0);
00055 filter1.setCriteria(ddv1,DDSpecificsFilter::equals);
00056 DDFilteredView fv1(cpv);
00057 fv1.addFilter(filter1);
00058 wcNames = getNames(fv1);
00059 edm::LogInfo("HcalTBSim") << "HcalTB06BeamSD:: Names to be tested for "
00060 << attribute << " = " << value << ": "
00061 << wcNames.size() << " paths";
00062 for (unsigned int i=0; i<wcNames.size(); i++)
00063 edm::LogInfo("HcalTBSim") << "HcalTB06BeamSD:: (" << i << ") "
00064 << wcNames[i];
00065
00066
00067 attribute = "ReadOutName";
00068 DDSpecificsFilter filter2;
00069 DDValue ddv2(attribute,name,0);
00070 filter2.setCriteria(ddv2,DDSpecificsFilter::equals);
00071 DDFilteredView fv2(cpv);
00072 fv2.addFilter(filter2);
00073 bool dodet = fv2.firstChild();
00074
00075 std::vector<G4String> matNames;
00076 std::vector<int> nocc;
00077 while (dodet) {
00078 const DDLogicalPart & log = fv2.logicalPart();
00079 matName = DDSplit(log.material().name()).first;
00080 bool notIn = true;
00081 for (unsigned int i=0; i<matNames.size(); i++)
00082 if (matName == matNames[i]) {notIn = false; nocc[i]++;}
00083 if (notIn) {
00084 matNames.push_back(matName);
00085 nocc.push_back(0);
00086 }
00087 dodet = fv2.next();
00088 }
00089 if (matNames.size() > 0) {
00090 matName = matNames[0];
00091 int occ = nocc[0];
00092 for (unsigned int i = 0; i < matNames.size(); i++) {
00093 if (nocc[i] > occ) {
00094 occ = nocc[i];
00095 matName = matNames[i];
00096 }
00097 }
00098 } else {
00099 matName = "Not Found";
00100 }
00101
00102 edm::LogInfo("HcalTBSim") << "HcalTB06BeamSD: Material name for "
00103 << attribute << " = " << name << ":" << matName;
00104 }
00105
00106 HcalTB06BeamSD::~HcalTB06BeamSD() {}
00107
00108 double HcalTB06BeamSD::getEnergyDeposit(G4Step* aStep) {
00109
00110 double destep = aStep->GetTotalEnergyDeposit();
00111 double weight = 1;
00112 if (useBirk) {
00113 G4Material* mat = aStep->GetPreStepPoint()->GetMaterial();
00114 if (mat->GetName() == matName)
00115 weight *= getAttenuation(aStep, birk1, birk2, birk3);
00116 }
00117 LogDebug("HcalTBSim") << "HcalTB06BeamSD: Detector "
00118 << aStep->GetPreStepPoint()->GetTouchable()->GetVolume()->GetName()
00119 << " weight " << weight;
00120 return weight*destep;
00121 }
00122
00123 uint32_t HcalTB06BeamSD::setDetUnitId(G4Step * aStep) {
00124
00125 G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
00126 const G4VTouchable* touch = preStepPoint->GetTouchable();
00127 G4String name = preStepPoint->GetPhysicalVolume()->GetName();
00128
00129 int det = 1;
00130 int lay = 0, x = 0, y = 0;
00131 if (!isItWireChamber(name)) {
00132 lay = (touch->GetReplicaNumber(0));
00133 } else {
00134 det = 2;
00135 lay = (touch->GetReplicaNumber(1));
00136 G4ThreeVector hitPoint = preStepPoint->GetPosition();
00137 G4ThreeVector localPoint = setToLocal(hitPoint, touch);
00138 x = (int)(localPoint.x()/(0.2*mm));
00139 y = (int)(localPoint.y()/(0.2*mm));
00140 }
00141
00142 return packIndex (det, lay, x, y);
00143 }
00144
00145 uint32_t HcalTB06BeamSD::packIndex(int det, int lay, int x, int y) {
00146
00147 int ix = 0, ixx = x;
00148 if (x < 0) { ix = 1; ixx =-x;}
00149 int iy = 0, iyy = y;
00150 if (y < 0) { iy = 1; iyy =-y;}
00151 uint32_t idx = (det&15)<<28;
00152 idx += (lay&127)<<21;
00153 idx += (iy&1)<<19;
00154 idx += (iyy&511)<<10;
00155 idx += (ix&1)<<9;
00156 idx += (ixx&511);
00157
00158 LogDebug("HcalTBSim") << "HcalTB06BeamSD: Detector " << det << " Layer "
00159 << lay << " x " << x << " " << ix << " " << ixx
00160 << " y " << y << " " << iy << " " << iyy << " ID "
00161 << std::hex << idx << std::dec;
00162 return idx;
00163 }
00164
00165 void HcalTB06BeamSD::unpackIndex(const uint32_t & idx, int& det, int& lay,
00166 int& x, int& y) {
00167
00168 det = (idx>>28)&15;
00169 lay = (idx>>21)&127;
00170 y = (idx>>10)&511; if (((idx>>19)&1) == 1) y = -y;
00171 x = (idx)&511; if (((idx>>9)&1) == 1) x = -x;
00172
00173 }
00174
00175 std::vector<G4String> HcalTB06BeamSD::getNames(DDFilteredView& fv) {
00176
00177 std::vector<G4String> tmp;
00178 bool dodet = fv.firstChild();
00179 while (dodet) {
00180 const DDLogicalPart & log = fv.logicalPart();
00181 G4String name = DDSplit(log.name()).first;
00182 bool ok = true;
00183 for (unsigned int i=0; i<tmp.size(); i++)
00184 if (name == tmp[i]) ok = false;
00185 if (ok) tmp.push_back(name);
00186 dodet = fv.next();
00187 }
00188 return tmp;
00189 }
00190
00191 bool HcalTB06BeamSD::isItWireChamber (G4String name) {
00192
00193 std::vector<G4String>::const_iterator it = wcNames.begin();
00194 for (; it != wcNames.end(); it++)
00195 if (name == *it) return true;
00196 return false;
00197 }