#include <SimG4CMS/HcalTestBeam/interface/HcalTB02SD.h>
Public Member Functions | |
virtual double | getEnergyDeposit (G4Step *) |
HcalTB02SD (G4String, const DDCompactView &, SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *) | |
virtual uint32_t | setDetUnitId (G4Step *step) |
void | setNumberingScheme (HcalTB02NumberingScheme *scheme) |
virtual | ~HcalTB02SD () |
Private Member Functions | |
double | crystalLength (G4String) |
double | curve_LY (G4String &, G4StepPoint *) |
void | initMap (G4String, const DDCompactView &) |
Private Attributes | |
double | birk1 |
double | birk2 |
double | birk3 |
std::map< G4String, double > | lengthMap |
HcalTB02NumberingScheme * | numberingScheme |
bool | useBirk |
bool | useWeight |
Description: Stores hits of Test Beam 2002 calorimeters
Usage: Activation is done using the XML file by choosing HcalTB02SensitiveDetector
Definition at line 30 of file HcalTB02SD.h.
HcalTB02SD::HcalTB02SD | ( | G4String | name, |
const DDCompactView & | cpv, | ||
SensitiveDetectorCatalog & | clg, | ||
edm::ParameterSet const & | p, | ||
const SimTrackManager * | manager | ||
) |
Definition at line 35 of file HcalTB02SD.cc.
References birk1, birk2, birk3, g, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), initMap(), LogDebug, setNumberingScheme(), useBirk, and useWeight.
: CaloSD(name, cpv, clg, p, manager), numberingScheme(0) { edm::ParameterSet m_SD = p.getParameter<edm::ParameterSet>("HcalTB02SD"); useBirk= m_SD.getUntrackedParameter<bool>("UseBirkLaw",false); birk1 = m_SD.getUntrackedParameter<double>("BirkC1",0.013)*(g/(MeV*cm2)); birk2 = m_SD.getUntrackedParameter<double>("BirkC2",0.0568); birk3 = m_SD.getUntrackedParameter<double>("BirkC3",1.75); useWeight= true; HcalTB02NumberingScheme* scheme=0; if (name == "EcalHitsEB") { scheme = dynamic_cast<HcalTB02NumberingScheme*>(new HcalTB02XtalNumberingScheme()); useBirk = false; } else if (name == "HcalHits") { scheme = dynamic_cast<HcalTB02NumberingScheme*>(new HcalTB02HcalNumberingScheme()); useWeight= false; } else {edm::LogWarning("HcalTBSim") << "HcalTB02SD: ReadoutName " << name << " not supported\n";} if (scheme) setNumberingScheme(scheme); LogDebug("HcalTBSim") << "***************************************************" << "\n" << "* *" << "\n" << "* Constructing a HcalTB02SD with name " << GetName() << "\n" << "* *" << "\n" << "***************************************************" ; edm::LogInfo("HcalTBSim") << "HcalTB02SD:: Use of Birks law is set to " << useBirk << " with three constants kB = " << birk1 << ", C1 = " << birk2 << ", C2 = " << birk3; if (useWeight) initMap(name,cpv); }
HcalTB02SD::~HcalTB02SD | ( | ) | [virtual] |
Definition at line 78 of file HcalTB02SD.cc.
References numberingScheme.
{ if (numberingScheme) delete numberingScheme; }
double HcalTB02SD::crystalLength | ( | G4String | name | ) | [private] |
double HcalTB02SD::curve_LY | ( | G4String & | nameVolume, |
G4StepPoint * | stepPoint | ||
) | [private] |
Definition at line 153 of file HcalTB02SD.cc.
References crystalLength(), LogDebug, CaloSD::setToLocal(), and CommonMethods::weight().
Referenced by getEnergyDeposit().
{ double weight = 1.; G4ThreeVector localPoint = setToLocal(stepPoint->GetPosition(), stepPoint->GetTouchable()); double crlength = crystalLength(nameVolume); double dapd = 0.5 * crlength - localPoint.z(); if (dapd >= -0.1 || dapd <= crlength+0.1) { if (dapd <= 100.) weight = 1.05 - dapd * 0.0005; } else { edm::LogWarning("HcalTBSim") << "HcalTB02SD: light coll curve : wrong " << "distance to APD " << dapd <<" crlength = " << crlength << " crystal name = " <<nameVolume << " z of localPoint = " << localPoint.z() << " take weight = " << weight; } LogDebug("HcalTBSim") << "HcalTB02SD, light coll curve : " << dapd << " crlength = " << crlength << " crystal name = " << nameVolume << " z of localPoint = " << localPoint.z() << " take weight = " << weight; return weight; }
double HcalTB02SD::getEnergyDeposit | ( | G4Step * | aStep | ) | [virtual] |
Reimplemented from CaloSD.
Definition at line 86 of file HcalTB02SD.cc.
References birk1, birk2, birk3, curve_LY(), CaloSD::getAttenuation(), LogDebug, NULL, CaloSD::preStepPoint, useBirk, useWeight, and CommonMethods::weight().
{ if (aStep == NULL) { return 0; } else { preStepPoint = aStep->GetPreStepPoint(); G4String nameVolume = preStepPoint->GetPhysicalVolume()->GetName(); // take into account light collection curve for crystals double weight = 1.; if (useWeight) weight *= curve_LY(nameVolume, preStepPoint); if (useBirk) weight *= getAttenuation(aStep, birk1, birk2, birk3); double edep = aStep->GetTotalEnergyDeposit() * weight; LogDebug("HcalTBSim") << "HcalTB02SD:: " << nameVolume <<" Light Collection Efficiency " << weight << " Weighted Energy Deposit " << edep/MeV << " MeV"; return edep; } }
void HcalTB02SD::initMap | ( | G4String | sd, |
const DDCompactView & | cpv | ||
) | [private] |
Definition at line 119 of file HcalTB02SD.cc.
References DDFilteredView::addFilter(), ddtrap, DDSpecificsFilter::equals, alcazmumu_cfi::filter, DDFilteredView::firstChild(), i, lengthMap, LogDebug, DDFilteredView::logicalPart(), SensitiveDetector::name, DDName::name(), DDBase< N, C >::name(), DDFilteredView::next(), DDSolid::parameters(), DDSpecificsFilter::setCriteria(), DDSolid::shape(), and DDLogicalPart::solid().
Referenced by HcalTB02SD().
{ G4String attribute = "ReadOutName"; DDSpecificsFilter filter; DDValue ddv(attribute,sd,0); filter.setCriteria(ddv,DDSpecificsFilter::equals); DDFilteredView fv(cpv); fv.addFilter(filter); fv.firstChild(); bool dodet=true; while (dodet) { const DDSolid & sol = fv.logicalPart().solid(); const std::vector<double> & paras = sol.parameters(); G4String name = sol.name().name(); LogDebug("HcalTBSim") << "HcalTB02SD::initMap (for " << sd << "): Solid " << name << " Shape " << sol.shape() << " Parameter 0 = " << paras[0]; if (sol.shape() == ddtrap) { double dz = 2*paras[0]; lengthMap.insert(std::pair<G4String,double>(name,dz)); } dodet = fv.next(); } LogDebug("HcalTBSim") << "HcalTB02SD: Length Table for " << attribute << " = " << sd << ":"; std::map<G4String,double>::const_iterator it = lengthMap.begin(); int i=0; for (; it != lengthMap.end(); it++, i++) { LogDebug("HcalTBSim") << " " << i << " " << it->first << " L = " << it->second; } }
uint32_t HcalTB02SD::setDetUnitId | ( | G4Step * | step | ) | [virtual] |
Implements CaloSD.
Definition at line 106 of file HcalTB02SD.cc.
References HcalTB02NumberingScheme::getUnitID(), and numberingScheme.
{ return (numberingScheme == 0 ? 0 : (uint32_t)(numberingScheme->getUnitID(aStep))); }
void HcalTB02SD::setNumberingScheme | ( | HcalTB02NumberingScheme * | scheme | ) |
Definition at line 110 of file HcalTB02SD.cc.
References numberingScheme.
Referenced by HcalTB02SD().
{ if (scheme != 0) { edm::LogInfo("HcalTBSim") << "HcalTB02SD: updates numbering scheme for " << GetName(); if (numberingScheme) delete numberingScheme; numberingScheme = scheme; } }
double HcalTB02SD::birk1 [private] |
Definition at line 51 of file HcalTB02SD.h.
Referenced by getEnergyDeposit(), and HcalTB02SD().
double HcalTB02SD::birk2 [private] |
Definition at line 51 of file HcalTB02SD.h.
Referenced by getEnergyDeposit(), and HcalTB02SD().
double HcalTB02SD::birk3 [private] |
Definition at line 51 of file HcalTB02SD.h.
Referenced by getEnergyDeposit(), and HcalTB02SD().
std::map<G4String,double> HcalTB02SD::lengthMap [private] |
Definition at line 52 of file HcalTB02SD.h.
Referenced by crystalLength(), and initMap().
Definition at line 48 of file HcalTB02SD.h.
Referenced by setDetUnitId(), setNumberingScheme(), and ~HcalTB02SD().
bool HcalTB02SD::useBirk [private] |
Definition at line 50 of file HcalTB02SD.h.
Referenced by getEnergyDeposit(), and HcalTB02SD().
bool HcalTB02SD::useWeight [private] |
Definition at line 49 of file HcalTB02SD.h.
Referenced by getEnergyDeposit(), and HcalTB02SD().