Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017 #include "SimG4CMS/HcalTestBeam/interface/HcalTB02SD.h"
00018 #include "SimG4CMS/HcalTestBeam/interface/HcalTB02HcalNumberingScheme.h"
00019 #include "SimG4CMS/HcalTestBeam/interface/HcalTB02XtalNumberingScheme.h"
00020 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00021 #include "DetectorDescription/Core/interface/DDFilter.h"
00022 #include "DetectorDescription/Core/interface/DDFilteredView.h"
00023 #include "DetectorDescription/Core/interface/DDSolid.h"
00024 #include "DetectorDescription/Core/interface/DDSplit.h"
00025 #include "DetectorDescription/Core/interface/DDValue.h"
00026
00027 #include "G4Step.hh"
00028 #include "G4Track.hh"
00029 #include "G4VProcess.hh"
00030
00031
00032
00033
00034
00035 HcalTB02SD::HcalTB02SD(G4String name, const DDCompactView & cpv,
00036 SensitiveDetectorCatalog & clg,
00037 edm::ParameterSet const & p,
00038 const SimTrackManager* manager) :
00039 CaloSD(name, cpv, clg, p, manager), numberingScheme(0) {
00040
00041 edm::ParameterSet m_SD = p.getParameter<edm::ParameterSet>("HcalTB02SD");
00042 useBirk= m_SD.getUntrackedParameter<bool>("UseBirkLaw",false);
00043 birk1 = m_SD.getUntrackedParameter<double>("BirkC1",0.013)*(g/(MeV*cm2));
00044 birk2 = m_SD.getUntrackedParameter<double>("BirkC2",0.0568);
00045 birk3 = m_SD.getUntrackedParameter<double>("BirkC3",1.75);
00046 useWeight= true;
00047
00048 HcalTB02NumberingScheme* scheme=0;
00049 if (name == "EcalHitsEB") {
00050 scheme = dynamic_cast<HcalTB02NumberingScheme*>(new HcalTB02XtalNumberingScheme());
00051 useBirk = false;
00052 } else if (name == "HcalHits") {
00053 scheme = dynamic_cast<HcalTB02NumberingScheme*>(new HcalTB02HcalNumberingScheme());
00054 useWeight= false;
00055 } else {edm::LogWarning("HcalTBSim") << "HcalTB02SD: ReadoutName " << name
00056 << " not supported\n";}
00057
00058 if (scheme) setNumberingScheme(scheme);
00059 LogDebug("HcalTBSim")
00060 << "***************************************************"
00061 << "\n"
00062 << "* *"
00063 << "\n"
00064 << "* Constructing a HcalTB02SD with name " << GetName()
00065 << "\n"
00066 << "* *"
00067 << "\n"
00068 << "***************************************************" ;
00069 edm::LogInfo("HcalTBSim") << "HcalTB02SD:: Use of Birks law is set to "
00070 << useBirk << " with three constants kB = "
00071 << birk1 << ", C1 = " << birk2 << ", C2 = "
00072 << birk3;
00073
00074 if (useWeight) initMap(name,cpv);
00075
00076 }
00077
00078 HcalTB02SD::~HcalTB02SD() {
00079 if (numberingScheme) delete numberingScheme;
00080 }
00081
00082
00083
00084
00085
00086 double HcalTB02SD::getEnergyDeposit(G4Step * aStep) {
00087
00088 if (aStep == NULL) {
00089 return 0;
00090 } else {
00091 preStepPoint = aStep->GetPreStepPoint();
00092 G4String nameVolume = preStepPoint->GetPhysicalVolume()->GetName();
00093
00094
00095 double weight = 1.;
00096 if (useWeight) weight *= curve_LY(nameVolume, preStepPoint);
00097 if (useBirk) weight *= getAttenuation(aStep, birk1, birk2, birk3);
00098 double edep = aStep->GetTotalEnergyDeposit() * weight;
00099 LogDebug("HcalTBSim") << "HcalTB02SD:: " << nameVolume
00100 <<" Light Collection Efficiency " << weight
00101 << " Weighted Energy Deposit " << edep/MeV << " MeV";
00102 return edep;
00103 }
00104 }
00105
00106 uint32_t HcalTB02SD::setDetUnitId(G4Step * aStep) {
00107 return (numberingScheme == 0 ? 0 : (uint32_t)(numberingScheme->getUnitID(aStep)));
00108 }
00109
00110 void HcalTB02SD::setNumberingScheme(HcalTB02NumberingScheme* scheme) {
00111 if (scheme != 0) {
00112 edm::LogInfo("HcalTBSim") << "HcalTB02SD: updates numbering scheme for "
00113 << GetName();
00114 if (numberingScheme) delete numberingScheme;
00115 numberingScheme = scheme;
00116 }
00117 }
00118
00119 void HcalTB02SD::initMap(G4String sd, const DDCompactView & cpv) {
00120
00121 G4String attribute = "ReadOutName";
00122 DDSpecificsFilter filter;
00123 DDValue ddv(attribute,sd,0);
00124 filter.setCriteria(ddv,DDSpecificsFilter::equals);
00125 DDFilteredView fv(cpv);
00126 fv.addFilter(filter);
00127 fv.firstChild();
00128
00129 bool dodet=true;
00130 while (dodet) {
00131 const DDSolid & sol = fv.logicalPart().solid();
00132 const std::vector<double> & paras = sol.parameters();
00133 G4String name = sol.name().name();
00134 LogDebug("HcalTBSim") << "HcalTB02SD::initMap (for " << sd << "): Solid "
00135 << name << " Shape " << sol.shape()
00136 << " Parameter 0 = " << paras[0];
00137 if (sol.shape() == ddtrap) {
00138 double dz = 2*paras[0];
00139 lengthMap.insert(std::pair<G4String,double>(name,dz));
00140 }
00141 dodet = fv.next();
00142 }
00143 LogDebug("HcalTBSim") << "HcalTB02SD: Length Table for " << attribute
00144 << " = " << sd << ":";
00145 std::map<G4String,double>::const_iterator it = lengthMap.begin();
00146 int i=0;
00147 for (; it != lengthMap.end(); it++, i++) {
00148 LogDebug("HcalTBSim") << " " << i << " " << it->first << " L = "
00149 << it->second;
00150 }
00151 }
00152
00153 double HcalTB02SD::curve_LY(G4String& nameVolume, G4StepPoint* stepPoint) {
00154
00155 double weight = 1.;
00156 G4ThreeVector localPoint = setToLocal(stepPoint->GetPosition(),
00157 stepPoint->GetTouchable());
00158 double crlength = crystalLength(nameVolume);
00159 double dapd = 0.5 * crlength - localPoint.z();
00160 if (dapd >= -0.1 || dapd <= crlength+0.1) {
00161 if (dapd <= 100.)
00162 weight = 1.05 - dapd * 0.0005;
00163 } else {
00164 edm::LogWarning("HcalTBSim") << "HcalTB02SD: light coll curve : wrong "
00165 << "distance to APD " << dapd <<" crlength = "
00166 << crlength << " crystal name = " <<nameVolume
00167 << " z of localPoint = " << localPoint.z()
00168 << " take weight = " << weight;
00169 }
00170 LogDebug("HcalTBSim") << "HcalTB02SD, light coll curve : " << dapd
00171 << " crlength = " << crlength
00172 << " crystal name = " << nameVolume
00173 << " z of localPoint = " << localPoint.z()
00174 << " take weight = " << weight;
00175 return weight;
00176 }
00177
00178 double HcalTB02SD::crystalLength(G4String name) {
00179
00180 double length = 230.;
00181 std::map<G4String,double>::const_iterator it = lengthMap.find(name);
00182 if (it != lengthMap.end()) length = it->second;
00183 return length;
00184 }