#include <SimG4CMS/Calo/interface/HCalSD.h>
Public Member Functions | |
virtual double | getEnergyDeposit (G4Step *) |
HCalSD (G4String, const DDCompactView &, SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *) | |
virtual bool | ProcessHits (G4Step *, G4TouchableHistory *) |
virtual uint32_t | setDetUnitId (G4Step *step) |
void | setNumberingScheme (HcalNumberingScheme *) |
virtual | ~HCalSD () |
Protected Member Functions | |
virtual void | initRun () |
Private Member Functions | |
std::vector< double > | getDDDArray (const std::string &, const DDsvalues_type &) |
void | getFromLibrary (G4Step *step) |
void | getFromParam (G4Step *step) |
void | getHitPMT (G4Step *step) |
std::vector< G4String > | getNames (DDFilteredView &) |
void | hitForFibre (G4Step *step) |
bool | isItFibre (G4String) |
bool | isItFibre (G4LogicalVolume *) |
bool | isItHF (G4String) |
bool | isItHF (G4Step *) |
bool | isItPMT (G4LogicalVolume *) |
bool | isItScintillator (G4Material *) |
double | layerWeight (int, G4ThreeVector, int, int) |
void | readWeightFromFile (std::string) |
uint32_t | setDetUnitId (int, G4ThreeVector, int, int) |
int | setTrackID (G4Step *step) |
Private Attributes | |
double | betaThr |
double | birk1 |
double | birk2 |
double | birk3 |
std::vector< G4LogicalVolume * > | fibreLV |
std::vector< G4String > | fibreNames |
std::vector< int > | hfLevels |
std::vector< G4LogicalVolume * > | hfLV |
std::vector< G4String > | hfNames |
HFShower * | hfshower |
std::vector< double > | layer0wt |
std::map< uint32_t, double > | layerWeights |
std::vector< G4Material * > | materials |
std::vector< G4String > | matNames |
G4int | mumPDG |
G4int | mupPDG |
HcalNumberingFromDDD * | numberingFromDDD |
HcalNumberingScheme * | numberingScheme |
std::vector< G4LogicalVolume * > | pmtLV |
std::vector< G4String > | pmtNames |
HFShowerLibrary * | showerLibrary |
HFShowerParam * | showerParam |
HFShowerPMT * | showerPMT |
bool | useBirk |
bool | useHF |
bool | useLayerWt |
bool | useParam |
bool | usePMTHit |
bool | useShowerLibrary |
Definition at line 29 of file HCalSD.h.
HCalSD::HCalSD | ( | G4String | name, | |
const DDCompactView & | cpv, | |||
SensitiveDetectorCatalog & | clg, | |||
edm::ParameterSet const & | p, | |||
const SimTrackManager * | manager | |||
) |
Definition at line 28 of file HCalSD.cc.
References DDFilteredView::addFilter(), betaThr, birk1, birk2, birk3, DDSplit(), DDSpecificsFilter::equals, fibreLV, fibreNames, file, first, DDFilteredView::firstChild(), g, getDDDArray(), SensitiveDetector::getNames(), edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), hfLevels, hfLV, hfNames, hfshower, i, isItFibre(), isItHF(), it, CaloSD::kmaxIon, CaloSD::kmaxNeutron, CaloSD::kmaxProton, layer0wt, level, funct::log(), LogDebug, DDFilteredView::logicalPart(), lv, DDLogicalPart::material(), materials, matNames, DDFilteredView::mergedSpecifics(), mumPDG, mupPDG, DDBase< N, C >::name(), DDFilteredView::next(), numberingFromDDD, pmtLV, pmtNames, readWeightFromFile(), DDSpecificsFilter::setCriteria(), setNumberingScheme(), showerLibrary, showerParam, showerPMT, CaloSD::suppressHeavy, sv, pyDBSRunClass::temp, useBirk, useHF, useLayerWt, useParam, usePMTHit, useShowerLibrary, and value.
00030 : 00031 CaloSD(name, cpv, clg, p, manager), numberingFromDDD(0), numberingScheme(0), 00032 showerLibrary(0), hfshower(0), showerParam(0), showerPMT(0) { 00033 00034 //static SimpleConfigurable<bool> on1(false, "HCalSD:UseBirkLaw"); 00035 //static SimpleConfigurable<double> bk1(0.013, "HCalSD:BirkC1"); 00036 //static SimpleConfigurable<double> bk2(0.0568,"HCalSD:BirkC2"); 00037 //static SimpleConfigurable<double> bk3(1.75, "HCalSD:BirkC3"); 00038 // Values from NIM 80 (1970) 239-244: as implemented in Geant3 00039 //static SimpleConfigurable<bool> on2(true,"HCalSD:UseShowerLibrary"); 00040 00041 edm::ParameterSet m_HC = p.getParameter<edm::ParameterSet>("HCalSD"); 00042 useBirk = m_HC.getParameter<bool>("UseBirkLaw"); 00043 birk1 = m_HC.getParameter<double>("BirkC1")*(g/(MeV*cm2)); 00044 birk2 = m_HC.getParameter<double>("BirkC2"); 00045 birk3 = m_HC.getParameter<double>("BirkC3"); 00046 useShowerLibrary = m_HC.getParameter<bool>("UseShowerLibrary"); 00047 useParam = m_HC.getParameter<bool>("UseParametrize"); 00048 bool testNumber = m_HC.getParameter<bool>("TestNumberingScheme"); 00049 usePMTHit = m_HC.getParameter<bool>("UsePMTHits"); 00050 betaThr = m_HC.getParameter<double>("BetaThreshold"); 00051 useHF = m_HC.getUntrackedParameter<bool>("UseHF",true); 00052 bool forTBH2 = m_HC.getUntrackedParameter<bool>("ForTBH2",false); 00053 useLayerWt = m_HC.getUntrackedParameter<bool>("UseLayerWt",false); 00054 std::string file = m_HC.getUntrackedParameter<std::string>("WtFile","None"); 00055 00056 LogDebug("HcalSim") << "***************************************************" 00057 << "\n" 00058 << "* *" 00059 << "\n" 00060 << "* Constructing a HCalSD with name " << name << "\n" 00061 << "* *" 00062 << "\n" 00063 << "***************************************************"; 00064 00065 edm::LogInfo("HcalSim") << "HCalSD:: Use of HF code is set to " << useHF 00066 << "\nUse of shower parametrization set to " 00067 << useParam << "\nUse of shower library is set to " 00068 << useShowerLibrary << "\nUse PMT Hit is set to " 00069 << usePMTHit << " with beta Threshold "<< betaThr 00070 << "\n Use of Birks law is set to " 00071 << useBirk << " with three constants kB = " 00072 << birk1 << ", C1 = " << birk2 << ", C2 = " << birk3; 00073 edm::LogInfo("HcalSim") << "HCalSD:: Suppression Flag " << suppressHeavy 00074 << " protons below " << kmaxProton << " MeV," 00075 << " neutrons below " << kmaxNeutron << " MeV and" 00076 << " ions below " << kmaxIon << " MeV"; 00077 00078 numberingFromDDD = new HcalNumberingFromDDD(name, cpv); 00079 HcalNumberingScheme* scheme; 00080 if (testNumber || forTBH2) 00081 scheme = dynamic_cast<HcalNumberingScheme*>(new HcalTestNumberingScheme(forTBH2)); 00082 else 00083 scheme = new HcalNumberingScheme(); 00084 setNumberingScheme(scheme); 00085 00086 const G4LogicalVolumeStore * lvs = G4LogicalVolumeStore::GetInstance(); 00087 std::vector<G4LogicalVolume *>::const_iterator lvcite; 00088 G4LogicalVolume* lv; 00089 std::string attribute, value; 00090 if (useHF) { 00091 if (useParam) showerParam = new HFShowerParam(name, cpv, p); 00092 else { 00093 if (useShowerLibrary) showerLibrary = new HFShowerLibrary(name, cpv, p); 00094 hfshower = new HFShower(name, cpv, p); 00095 } 00096 00097 // HF volume names 00098 attribute = "Volume"; 00099 value = "HF"; 00100 DDSpecificsFilter filter0; 00101 DDValue ddv0(attribute,value,0); 00102 filter0.setCriteria(ddv0,DDSpecificsFilter::equals); 00103 DDFilteredView fv0(cpv); 00104 fv0.addFilter(filter0); 00105 hfNames = getNames(fv0); 00106 fv0.firstChild(); 00107 DDsvalues_type sv0(fv0.mergedSpecifics()); 00108 std::vector<double> temp = getDDDArray("Levels",sv0); 00109 edm::LogInfo("HcalSim") << "HCalSD: Names to be tested for " << attribute 00110 << " = " << value << " has " << hfNames.size() 00111 << " elements"; 00112 for (unsigned int i=0; i<hfNames.size(); i++) { 00113 G4String namv = hfNames[i]; 00114 lv = 0; 00115 for (lvcite = lvs->begin(); lvcite != lvs->end(); lvcite++) 00116 if ((*lvcite)->GetName() == namv) { 00117 lv = (*lvcite); 00118 break; 00119 } 00120 hfLV.push_back(lv); 00121 int level = static_cast<int>(temp[i]); 00122 hfLevels.push_back(level); 00123 edm::LogInfo("HcalSim") << "HCalSD: HF[" << i << "] = " << hfNames[i] 00124 << " LV " << hfLV[i] << " at level " 00125 << hfLevels[i]; 00126 } 00127 00128 // HF Fibre volume names 00129 value = "HFFibre"; 00130 DDSpecificsFilter filter1; 00131 DDValue ddv1(attribute,value,0); 00132 filter1.setCriteria(ddv1,DDSpecificsFilter::equals); 00133 DDFilteredView fv1(cpv); 00134 fv1.addFilter(filter1); 00135 fibreNames = getNames(fv1); 00136 edm::LogInfo("HcalSim") << "HCalSD: Names to be tested for " << attribute 00137 << " = " << value << ":"; 00138 for (unsigned int i=0; i<fibreNames.size(); i++) { 00139 G4String namv = fibreNames[i]; 00140 lv = 0; 00141 for (lvcite = lvs->begin(); lvcite != lvs->end(); lvcite++) 00142 if ((*lvcite)->GetName() == namv) { 00143 lv = (*lvcite); 00144 break; 00145 } 00146 fibreLV.push_back(lv); 00147 edm::LogInfo("HcalSim") << "HCalSD: (" << i << ") " << fibreNames[i] 00148 << " LV " << fibreLV[i]; 00149 } 00150 00151 // HF PMT volume names 00152 value = "HFPMT"; 00153 DDSpecificsFilter filter3; 00154 DDValue ddv3(attribute,value,0); 00155 filter3.setCriteria(ddv3,DDSpecificsFilter::equals); 00156 DDFilteredView fv3(cpv); 00157 fv3.addFilter(filter3); 00158 pmtNames = getNames(fv3); 00159 edm::LogInfo("HcalSim") << "HCalSD: Names to be tested for " << attribute 00160 << " = " << value << " have " << pmtNames.size() 00161 << " entries"; 00162 for (unsigned int i=0; i<pmtNames.size(); i++) { 00163 G4String namv = pmtNames[i]; 00164 lv = 0; 00165 for (lvcite = lvs->begin(); lvcite != lvs->end(); lvcite++) 00166 if ((*lvcite)->GetName() == namv) { 00167 lv = (*lvcite); 00168 break; 00169 } 00170 pmtLV.push_back(lv); 00171 edm::LogInfo("HcalSim") << "HCalSD: (" << i << ") " << pmtNames[i] 00172 << " LV " << pmtLV[i]; 00173 } 00174 if (pmtNames.size() > 0) showerPMT = new HFShowerPMT (name, cpv, p); 00175 } 00176 00177 //Material list for HB/HE/HO sensitive detectors 00178 attribute = "ReadOutName"; 00179 DDSpecificsFilter filter2; 00180 DDValue ddv2(attribute,name,0); 00181 filter2.setCriteria(ddv2,DDSpecificsFilter::equals); 00182 DDFilteredView fv2(cpv); 00183 fv2.addFilter(filter2); 00184 bool dodet = fv2.firstChild(); 00185 00186 DDsvalues_type sv(fv2.mergedSpecifics()); 00187 //Layer0 Weight 00188 layer0wt = getDDDArray("Layer0Wt",sv); 00189 edm::LogInfo("HcalSim") << "HCalSD: " << layer0wt.size() << " Layer0Wt"; 00190 for (unsigned int it=0; it<layer0wt.size(); it++) 00191 edm::LogInfo("HcalSim") << "HCalSD: [" << it << "] " << layer0wt[it]; 00192 00193 const G4MaterialTable * matTab = G4Material::GetMaterialTable(); 00194 std::vector<G4Material*>::const_iterator matite; 00195 while (dodet) { 00196 const DDLogicalPart & log = fv2.logicalPart(); 00197 G4String namx = DDSplit(log.name()).first; 00198 if (!isItHF(namx) && !isItFibre(namx)) { 00199 namx = DDSplit(log.material().name()).first; 00200 bool notIn = true; 00201 for (unsigned int i=0; i<matNames.size(); i++) 00202 if (namx == matNames[i]) notIn = false; 00203 if (notIn) { 00204 matNames.push_back(namx); 00205 G4Material* mat; 00206 for (matite = matTab->begin(); matite != matTab->end(); matite++) 00207 if ((*matite)->GetName() == namx) { 00208 mat = (*matite); 00209 break; 00210 } 00211 materials.push_back(mat); 00212 } 00213 } 00214 dodet = fv2.next(); 00215 } 00216 00217 edm::LogInfo("HcalSim") << "HCalSD: Material names for " << attribute 00218 << " = " << name << ":"; 00219 for (unsigned int i=0; i<matNames.size(); i++) 00220 edm::LogInfo("HcalSim") << "HCalSD: (" << i << ") " << matNames[i] 00221 << " pointer " << materials[i]; 00222 00223 mumPDG = mupPDG = 0; 00224 00225 if (useLayerWt) readWeightFromFile(file); 00226 }
HCalSD::~HCalSD | ( | ) | [virtual] |
Definition at line 228 of file HCalSD.cc.
References hfshower, numberingFromDDD, numberingScheme, showerLibrary, showerParam, and showerPMT.
00228 { 00229 if (numberingFromDDD) delete numberingFromDDD; 00230 if (numberingScheme) delete numberingScheme; 00231 if (showerLibrary) delete showerLibrary; 00232 if (hfshower) delete hfshower; 00233 if (showerParam) delete showerParam; 00234 if (showerPMT) delete showerPMT; 00235 }
std::vector< double > HCalSD::getDDDArray | ( | const std::string & | str, | |
const DDsvalues_type & | sv | |||
) | [private] |
Definition at line 386 of file HCalSD.cc.
References DDfetch(), DDValue::doubles(), Exception, LogDebug, and value.
Referenced by HCalSD().
00387 { 00388 00389 LogDebug("HcalSim") << "HCalSD:getDDDArray called for " << str; 00390 DDValue value(str); 00391 if (DDfetch(&sv,value)) { 00392 LogDebug("HcalSim") << value; 00393 const std::vector<double> & fvec = value.doubles(); 00394 int nval = fvec.size(); 00395 if (nval < 1) { 00396 edm::LogError("HcalSim") << "HCalSD : # of " << str << " bins " << nval 00397 << " < 2 ==> illegal"; 00398 throw cms::Exception("Unknown", "HCalSD") 00399 << "nval < 2 for array " << str << "\n"; 00400 } 00401 00402 return fvec; 00403 } else { 00404 edm::LogError("HcalSim") << "HCalSD : cannot get array " << str; 00405 throw cms::Exception("Unknown", "HCalSD") 00406 << "cannot get array " << str << "\n"; 00407 } 00408 }
double HCalSD::getEnergyDeposit | ( | G4Step * | aStep | ) | [virtual] |
Reimplemented from CaloSD.
Definition at line 293 of file HCalSD.cc.
References birk1, birk2, birk3, CaloSD::getAttenuation(), isItScintillator(), TrackInformation::isPrimary(), CaloSD::kmaxIon, CaloSD::kmaxNeutron, CaloSD::kmaxProton, layer0wt, layerWeight(), LogDebug, CaloSD::suppressHeavy, CaloSD::theTrack, useBirk, useLayerWt, and weight.
00293 { 00294 00295 double destep = aStep->GetTotalEnergyDeposit(); 00296 double weight = 1; 00297 const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable(); 00298 int depth = (touch->GetReplicaNumber(0))%10; 00299 int det = ((touch->GetReplicaNumber(1))/1000)-3; 00300 if (depth==0 && (det==0 || det==1)) weight = layer0wt[det]; 00301 if (useLayerWt) { 00302 int lay = (touch->GetReplicaNumber(0)/10)%100 + 1; 00303 G4ThreeVector hitPoint = aStep->GetPreStepPoint()->GetPosition(); 00304 weight = layerWeight(det+3, hitPoint, depth, lay); 00305 } 00306 00307 if (suppressHeavy) { 00308 G4Track* theTrack = aStep->GetTrack(); 00309 TrackInformation * trkInfo = (TrackInformation *)(theTrack->GetUserInformation()); 00310 if (trkInfo) { 00311 int pdg = theTrack->GetDefinition()->GetPDGEncoding(); 00312 if (!(trkInfo->isPrimary())) { // Only secondary particles 00313 double ke = theTrack->GetKineticEnergy()/MeV; 00314 if (((pdg/1000000000 == 1 && ((pdg/10000)%100) > 0 && 00315 ((pdg/10)%100) > 0)) && (ke<kmaxIon)) weight = 0; 00316 if ((pdg == 2212) && (ke < kmaxProton)) weight = 0; 00317 if ((pdg == 2112) && (ke < kmaxNeutron)) weight = 0; 00318 if (weight == 0) { 00319 LogDebug("HcalSim") << "Ignore Track " << theTrack->GetTrackID() 00320 << " Type " << theTrack->GetDefinition()->GetParticleName() 00321 << " Kinetic Energy " << ke << " MeV"; 00322 } 00323 } 00324 } 00325 } 00326 00327 double weight0 = weight; 00328 if (useBirk) { 00329 G4Material* mat = aStep->GetPreStepPoint()->GetMaterial(); 00330 if (isItScintillator(mat)) 00331 weight *= getAttenuation(aStep, birk1, birk2, birk3); 00332 } 00333 LogDebug("HcalSim") << "HCalSD: Detector " << det+3 << " Depth " << depth 00334 << " weight " << weight0 << " " << weight; 00335 return weight*destep; 00336 }
void HCalSD::getFromLibrary | ( | G4Step * | step | ) | [private] |
Definition at line 479 of file HCalSD.cc.
References CaloSD::checkHit(), CaloSD::createNewHit(), CaloSD::currentHit, CaloSD::currentID, CaloSD::edepositEM, CaloSD::edepositHAD, CaloSD::emPDG, CaloSD::energyCut, CaloSD::epPDG, CaloSD::gammaPDG, HFShowerLibrary::getDepth(), HFShowerLibrary::getHits(), HFShowerLibrary::getPosHit(), HFShowerLibrary::getTSlice(), i, LogDebug, CaloSD::posGlobal, CaloSD::preStepPoint, CaloSD::previousID, CaloSD::resetForNewPrimary(), setDetUnitId(), CaloHitID::setID(), showerLibrary, CaloSD::theTrack, and CaloSD::updateHit().
Referenced by ProcessHits().
00479 { 00480 00481 preStepPoint = aStep->GetPreStepPoint(); 00482 theTrack = aStep->GetTrack(); 00483 00484 int nhit = showerLibrary->getHits(aStep); 00485 int det = 5; 00486 00487 double etrack = preStepPoint->GetKineticEnergy(); 00488 int primaryID = 0; 00489 if (etrack >= energyCut) { 00490 primaryID = theTrack->GetTrackID(); 00491 } else { 00492 primaryID = theTrack->GetParentID(); 00493 if (primaryID == 0) primaryID = theTrack->GetTrackID(); 00494 } 00495 00496 // Reset entry point for new primary 00497 posGlobal = preStepPoint->GetPosition(); 00498 resetForNewPrimary(posGlobal, etrack); 00499 // int primaryID = setTrackID(aStep); 00500 00501 G4int particleCode = theTrack->GetDefinition()->GetPDGEncoding(); 00502 if ( particleCode == emPDG || particleCode == epPDG || particleCode == gammaPDG ) { 00503 edepositEM = 1.*GeV; edepositHAD = 0.; 00504 } else { 00505 edepositEM = 0.; edepositHAD = 1.*GeV; 00506 } 00507 00508 LogDebug("HcalSim") << "HCalSD::getFromLibrary " << nhit << " hits for " 00509 << GetName() << " of " << primaryID << " with " 00510 << theTrack->GetDefinition()->GetParticleName() << " of " 00511 << preStepPoint->GetKineticEnergy()/GeV << " GeV"; 00512 00513 for (int i=0; i<nhit; i++) { 00514 G4ThreeVector hitPoint = showerLibrary->getPosHit(i); 00515 int depth = showerLibrary->getDepth(i); 00516 double time = showerLibrary->getTSlice(i); 00517 unsigned int unitID = setDetUnitId(det, hitPoint, depth); 00518 currentID.setID(unitID, time, primaryID, 0); 00519 00520 // check if it is in the same unit and timeslice as the previous one 00521 if (currentID == previousID) { 00522 updateHit(currentHit); 00523 } else { 00524 if (!checkHit()) currentHit = createNewHit(); 00525 } 00526 } 00527 00528 //Now kill the current track 00529 if (nhit >= 0) { 00530 theTrack->SetTrackStatus(fStopAndKill); 00531 } 00532 }
void HCalSD::getFromParam | ( | G4Step * | step | ) | [private] |
Definition at line 579 of file HCalSD.cc.
References CaloSD::checkHit(), CaloSD::createNewHit(), CaloSD::currentHit, CaloSD::currentID, CaloSD::edepositEM, CaloSD::edepositHAD, HFShowerParam::getDepth(), HFShowerParam::getHits(), HFShowerParam::getPosHit(), HFShowerParam::getTSlice(), i, LogDebug, CaloSD::posGlobal, CaloSD::preStepPoint, CaloSD::previousID, setDetUnitId(), CaloHitID::setID(), setTrackID(), showerParam, and CaloSD::updateHit().
Referenced by ProcessHits().
00579 { 00580 00581 std::vector<double> edeps = showerParam->getHits(aStep); 00582 int nHit = static_cast<int>(edeps.size()); 00583 00584 if (nHit > 0) { 00585 edepositEM = edeps[0]*GeV; 00586 edepositHAD = 0.; 00587 00588 preStepPoint = aStep->GetPreStepPoint(); 00589 int primaryID = setTrackID(aStep); 00590 00591 int det = 5; 00592 LogDebug("HcalSim") << "HCalSD::getFromParam " << nHit << " hits for " 00593 << GetName() << " of " << primaryID << " with " 00594 << aStep->GetTrack()->GetDefinition()->GetParticleName() 00595 << " of " << preStepPoint->GetKineticEnergy()/GeV 00596 << " GeV in detector type " << det; 00597 00598 for (int i = 0; i<nHit; i++) { 00599 G4ThreeVector hitPoint = showerParam->getPosHit(i); 00600 int depth = showerParam->getDepth(i); 00601 double time = showerParam->getTSlice(i); 00602 unsigned int unitID = setDetUnitId(det, hitPoint, depth); 00603 currentID.setID(unitID, time, primaryID, 0); 00604 00605 // check if it is in the same unit and timeslice as the previous one 00606 if (currentID == previousID) { 00607 updateHit(currentHit); 00608 } else { 00609 posGlobal = preStepPoint->GetPosition(); 00610 if (!checkHit()) currentHit = createNewHit(); 00611 } 00612 } 00613 } 00614 }
void HCalSD::getHitPMT | ( | G4Step * | step | ) | [private] |
Definition at line 616 of file HCalSD.cc.
References DeDxTools::beta(), betaThr, CaloSD::checkHit(), CaloSD::createNewHit(), CaloSD::currentHit, CaloSD::currentID, CaloSD::edepositEM, CaloSD::edepositHAD, CaloSD::energyCut, HFShowerPMT::getHits(), HFShowerPMT::getRadius(), HcalNumberingScheme::getUnitID(), HcalForward, LogDebug, numberingFromDDD, numberingScheme, phi, CaloSD::posGlobal, CaloSD::preStepPoint, CaloSD::previousID, CaloSD::resetForNewPrimary(), CaloHitID::setID(), showerPMT, CaloSD::theTrack, tmp, HcalNumberingFromDDD::unitID(), and CaloSD::updateHit().
Referenced by ProcessHits().
00616 { 00617 00618 preStepPoint = aStep->GetPreStepPoint(); 00619 theTrack = aStep->GetTrack(); 00620 double edep = showerPMT->getHits(aStep); 00621 00622 if (edep > 0) { 00623 double etrack = preStepPoint->GetKineticEnergy(); 00624 int primaryID = 0; 00625 if (etrack >= energyCut) { 00626 primaryID = theTrack->GetTrackID(); 00627 } else { 00628 primaryID = theTrack->GetParentID(); 00629 if (primaryID == 0) primaryID = theTrack->GetTrackID(); 00630 } 00631 // Reset entry point for new primary 00632 posGlobal = preStepPoint->GetPosition(); 00633 resetForNewPrimary(posGlobal, etrack); 00634 00635 // 00636 int det = static_cast<int>(HcalForward); 00637 G4ThreeVector hitPoint = preStepPoint->GetPosition(); 00638 double rr = (hitPoint.x()*hitPoint.x() + hitPoint.y()*hitPoint.y()); 00639 double phi = (rr == 0. ? 0. :atan2(hitPoint.y(),hitPoint.x())); 00640 double etaR = showerPMT->getRadius(); 00641 int depth = 1; 00642 if (etaR < 0) { 00643 depth = 2; 00644 etaR =-etaR; 00645 } 00646 if (hitPoint.z() < 0) etaR =-etaR; 00647 LogDebug("HcalSim") << "HCalSD::Hit for Detector " << det << " etaR " 00648 << etaR << " phi " << phi/deg << " depth " << depth; 00649 00650 double time = (aStep->GetPostStepPoint()->GetGlobalTime()); 00651 uint32_t unitID = 0; 00652 if (numberingFromDDD) { 00653 HcalNumberingFromDDD::HcalID tmp = numberingFromDDD->unitID(det,etaR,phi, 00654 depth,1); 00655 if (numberingScheme) unitID = numberingScheme->getUnitID(tmp); 00656 } 00657 currentID.setID(unitID, time, primaryID, 1); 00658 00659 double beta = preStepPoint->GetBeta(); 00660 if (beta > betaThr) { 00661 edepositEM = edep*GeV; edepositHAD = aStep->GetTotalEnergyDeposit(); 00662 } else { 00663 edepositEM = 0.; edepositHAD = aStep->GetTotalEnergyDeposit(); 00664 } 00665 LogDebug("HcalSim") << "HCalSD::getHitPMT 1 hit for " << GetName() 00666 << " of " << primaryID << " with " 00667 << theTrack->GetDefinition()->GetParticleName() 00668 << " of " << preStepPoint->GetKineticEnergy()/GeV 00669 << " GeV with velocity " << beta << " UnitID " 00670 << std::hex << unitID << std::dec; 00671 00672 // check if it is in the same unit and timeslice as the previous one 00673 if (currentID == previousID) { 00674 updateHit(currentHit); 00675 } else { 00676 if (!checkHit()) currentHit = createNewHit(); 00677 } 00678 } 00679 }
std::vector< G4String > HCalSD::getNames | ( | DDFilteredView & | fv | ) | [private] |
Definition at line 410 of file HCalSD.cc.
References DDSplit(), first, DDFilteredView::firstChild(), i, funct::log(), DDFilteredView::logicalPart(), DDBase< N, C >::name(), DDFilteredView::next(), and tmp.
00410 { 00411 00412 std::vector<G4String> tmp; 00413 bool dodet = fv.firstChild(); 00414 while (dodet) { 00415 const DDLogicalPart & log = fv.logicalPart(); 00416 G4String namx = DDSplit(log.name()).first; 00417 bool ok = true; 00418 for (unsigned int i=0; i<tmp.size(); i++) 00419 if (namx == tmp[i]) ok = false; 00420 if (ok) tmp.push_back(namx); 00421 dodet = fv.next(); 00422 } 00423 return tmp; 00424 }
void HCalSD::hitForFibre | ( | G4Step * | step | ) | [private] |
Definition at line 534 of file HCalSD.cc.
References CaloSD::checkHit(), CaloSD::createNewHit(), CaloSD::currentHit, CaloSD::currentID, CaloSD::edepositEM, CaloSD::edepositHAD, CaloSD::emPDG, CaloSD::epPDG, CaloSD::gammaPDG, HFShower::getHits(), HFShower::getTSlice(), hfshower, i, LogDebug, CaloSD::posGlobal, CaloSD::preStepPoint, CaloSD::previousID, setDetUnitId(), CaloHitID::setID(), setTrackID(), CaloSD::theTrack, and CaloSD::updateHit().
Referenced by ProcessHits().
00534 { 00535 00536 preStepPoint = aStep->GetPreStepPoint(); 00537 theTrack = aStep->GetTrack(); 00538 int primaryID = setTrackID(aStep); 00539 00540 int det = 5; 00541 int nHit = hfshower->getHits(aStep); 00542 00543 G4int particleCode = theTrack->GetDefinition()->GetPDGEncoding(); 00544 if ( particleCode == emPDG || particleCode == epPDG || particleCode == gammaPDG ) { 00545 edepositEM = 1.*GeV; edepositHAD = 0.; 00546 } else { 00547 edepositEM = 0.; edepositHAD = 1.*GeV; 00548 } 00549 00550 LogDebug("HcalSim") << "HCalSD::hitForFibre " << nHit << " hits for " 00551 << GetName() << " of " << primaryID << " with " 00552 << theTrack->GetDefinition()->GetParticleName() << " of " 00553 << preStepPoint->GetKineticEnergy()/GeV 00554 << " GeV in detector type " << det; 00555 00556 if (nHit > 0) { 00557 00558 G4ThreeVector hitPoint = preStepPoint->GetPosition(); 00559 int depth = 00560 (preStepPoint->GetTouchable()->GetReplicaNumber(0))%10; 00561 unsigned int unitID = setDetUnitId(det, hitPoint, depth); 00562 00563 for (int i=0; i<nHit; i++) { 00564 double time = hfshower->getTSlice(i); 00565 currentID.setID(unitID, time, primaryID, 0); 00566 00567 // check if it is in the same unit and timeslice as the previous one 00568 if (currentID == previousID) { 00569 updateHit(currentHit); 00570 } else { 00571 posGlobal = preStepPoint->GetPosition(); 00572 if (!checkHit()) currentHit = createNewHit(); 00573 } 00574 } 00575 } 00576 00577 }
void HCalSD::initRun | ( | ) | [protected, virtual] |
Reimplemented from CaloSD.
Definition at line 360 of file HCalSD.cc.
References HFShowerLibrary::initRun(), LogDebug, mumPDG, mupPDG, and showerLibrary.
00360 { 00361 00362 G4ParticleTable * theParticleTable = G4ParticleTable::GetParticleTable(); 00363 G4String particleName; 00364 mumPDG = theParticleTable->FindParticle(particleName="mu-")->GetPDGEncoding(); 00365 mupPDG = theParticleTable->FindParticle(particleName="mu+")->GetPDGEncoding(); 00366 LogDebug("HcalSim") << "HCalSD: Particle code for mu- = " << mumPDG 00367 << " for mu+ = " << mupPDG; 00368 00369 if (showerLibrary) showerLibrary->initRun(theParticleTable); 00370 }
bool HCalSD::isItFibre | ( | G4String | name | ) | [private] |
Definition at line 455 of file HCalSD.cc.
References fibreNames, and it.
00455 { 00456 00457 std::vector<G4String>::const_iterator it = fibreNames.begin(); 00458 for (; it != fibreNames.end(); it++) 00459 if (name == *it) return true; 00460 return false; 00461 }
bool HCalSD::isItFibre | ( | G4LogicalVolume * | lv | ) | [private] |
bool HCalSD::isItHF | ( | G4String | name | ) | [private] |
bool HCalSD::isItHF | ( | G4Step * | aStep | ) | [private] |
Definition at line 426 of file HCalSD.cc.
References hfLevels, hfLV, hfNames, it, and lv.
Referenced by HCalSD(), and ProcessHits().
00426 { 00427 00428 const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable(); 00429 int levels = ((touch->GetHistoryDepth())+1); 00430 for (unsigned int it=0; it < hfNames.size(); it++) { 00431 if (levels >= hfLevels[it]) { 00432 G4LogicalVolume* lv = touch->GetVolume(levels-hfLevels[it])->GetLogicalVolume(); 00433 if (lv == hfLV[it]) return true; 00434 } 00435 } 00436 return false; 00437 }
bool HCalSD::isItPMT | ( | G4LogicalVolume * | lv | ) | [private] |
bool HCalSD::isItScintillator | ( | G4Material * | mat | ) | [private] |
Definition at line 724 of file HCalSD.cc.
References HcalNumberingFromDDD::HcalID::etaR, HcalNumberingFromDDD::HcalID::lay, layerWeights, LogDebug, numberingFromDDD, HcalTestNumbering::packHcalIndex(), HcalNumberingFromDDD::HcalID::phis, HcalNumberingFromDDD::HcalID::subdet, tmp, HcalNumberingFromDDD::unitID(), and HcalNumberingFromDDD::HcalID::zside.
Referenced by getEnergyDeposit().
00724 { 00725 00726 double wt = 1.; 00727 if (numberingFromDDD) { 00728 //get the ID's as eta, phi, depth, ... indices 00729 HcalNumberingFromDDD::HcalID tmp = numberingFromDDD->unitID(det, pos, 00730 depth, lay); 00731 uint32_t id = HcalTestNumbering::packHcalIndex(tmp.subdet, tmp.zside, 00732 1, tmp.etaR, tmp.phis, 00733 tmp.lay); 00734 std::map<uint32_t,double>::const_iterator ite = layerWeights.find(id); 00735 if (ite != layerWeights.end()) wt = ite->second; 00736 LogDebug("HcalSim") << "HCalSD::layerWeight: ID " << std::hex << id 00737 << std::dec << " (" << tmp.subdet << "/" << tmp.zside 00738 << "/1/" << tmp.etaR << "/" << tmp.phis << "/" 00739 << tmp.lay << ") Weight " << wt; 00740 } 00741 return wt; 00742 }
bool HCalSD::ProcessHits | ( | G4Step * | aStep, | |
G4TouchableHistory * | ||||
) | [virtual] |
Reimplemented from CaloSD.
Definition at line 237 of file HCalSD.cc.
References CaloSD::createNewHit(), CaloSD::currentHit, CaloSD::edepositEM, CaloSD::edepositHAD, getFromLibrary(), getFromParam(), getHitPMT(), CaloSD::getStepInfo(), CaloSD::hitExists(), hitForFibre(), isItFibre(), isItHF(), isItPMT(), LogDebug, lv, mumPDG, mupPDG, SensitiveDetector::NaNTrap(), NULL, showerPMT, useParam, usePMTHit, and useShowerLibrary.
00237 { 00238 // TimeMe t1( theHitTimer, false); 00239 00240 NaNTrap( aStep ) ; 00241 00242 if (aStep == NULL) { 00243 return true; 00244 } else { 00245 G4LogicalVolume* lv = 00246 aStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume(); 00247 G4String nameVolume = lv->GetName(); 00248 if (isItHF(aStep)) { 00249 G4int parCode =aStep->GetTrack()->GetDefinition()->GetPDGEncoding(); 00250 if (useParam) { 00251 LogDebug("HcalSim") << "HCalSD: Hit from parametrization in " 00252 << nameVolume << " for Track " 00253 << aStep->GetTrack()->GetTrackID() 00254 <<" (" << aStep->GetTrack()->GetDefinition()->GetParticleName() << ")"; 00255 getFromParam(aStep); 00256 } else { 00257 bool notaMuon = true; 00258 if (parCode == mupPDG || parCode == mumPDG ) notaMuon = false; 00259 if (useShowerLibrary && notaMuon) { 00260 LogDebug("HcalSim") << "HCalSD: Starts shower library from " 00261 << nameVolume << " for Track " 00262 << aStep->GetTrack()->GetTrackID() 00263 <<" (" << aStep->GetTrack()->GetDefinition()->GetParticleName() << ")"; 00264 getFromLibrary(aStep); 00265 } else if (isItFibre(lv)) { 00266 LogDebug("HcalSim") << "HCalSD: Hit at Fibre in " << nameVolume 00267 << " for Track " 00268 << aStep->GetTrack()->GetTrackID() 00269 <<" (" << aStep->GetTrack()->GetDefinition()->GetParticleName() << ")"; 00270 hitForFibre(aStep); 00271 } 00272 } 00273 } else if (isItPMT(lv)) { 00274 LogDebug("HcalSim") << "HCalSD: Hit from PMT parametrization from " 00275 << nameVolume << " for Track " 00276 << aStep->GetTrack()->GetTrackID() << " (" 00277 << aStep->GetTrack()->GetDefinition()->GetParticleName() << ")"; 00278 if (usePMTHit && showerPMT) getHitPMT(aStep); 00279 } else { 00280 LogDebug("HcalSim") << "HCalSD: Hit from standard path from " 00281 << nameVolume << " for Track " 00282 << aStep->GetTrack()->GetTrackID() << " (" 00283 << aStep->GetTrack()->GetDefinition()->GetParticleName() << ")"; 00284 if (getStepInfo(aStep)) { 00285 if (hitExists() == false && edepositEM+edepositHAD>0.) 00286 currentHit = createNewHit(); 00287 } 00288 } 00289 return true; 00290 } 00291 }
void HCalSD::readWeightFromFile | ( | std::string | fName | ) | [private] |
Definition at line 700 of file HCalSD.cc.
References in, TreeToEdges::infile, layerWeights, LogDebug, HcalTestNumbering::packHcalIndex(), phi, and useLayerWt.
Referenced by HCalSD().
00700 { 00701 00702 std::ifstream infile; 00703 int entry=0; 00704 infile.open(fName.c_str(), std::ios::in); 00705 if (infile) { 00706 int det, zside, etaR, phi, lay; 00707 double wt; 00708 while (infile >> det >> zside >> etaR >> phi >> lay >> wt) { 00709 uint32_t id = HcalTestNumbering::packHcalIndex(det,zside,1,etaR,phi,lay); 00710 layerWeights.insert(std::pair<uint32_t,double>(id,wt)); 00711 entry++; 00712 LogDebug("HcalSim") << "HCalSD::readWeightFromFile:Entry " << entry 00713 << " ID " << std::hex << id << std::dec << " (" 00714 << det << "/" << zside << "/1/" << etaR << "/" 00715 << phi << "/" << lay << ") Weight " << wt; 00716 } 00717 infile.close(); 00718 } 00719 edm::LogInfo("HcalSim") << "HCalSD::readWeightFromFile: reads " << entry 00720 << " weights from " << fName; 00721 if (entry <= 0) useLayerWt = false; 00722 }
Definition at line 372 of file HCalSD.cc.
References HcalNumberingScheme::getUnitID(), numberingFromDDD, numberingScheme, tmp, and HcalNumberingFromDDD::unitID().
00373 { 00374 00375 uint32_t id = 0; 00376 if (numberingFromDDD) { 00377 //get the ID's as eta, phi, depth, ... indices 00378 HcalNumberingFromDDD::HcalID tmp = numberingFromDDD->unitID(det, pos, 00379 depth, lay); 00380 //get the ID 00381 if (numberingScheme) id = numberingScheme->getUnitID(tmp); 00382 } 00383 return id; 00384 }
uint32_t HCalSD::setDetUnitId | ( | G4Step * | step | ) | [virtual] |
Implements CaloSD.
Definition at line 338 of file HCalSD.cc.
References CaloSD::preStepPoint.
Referenced by getFromLibrary(), getFromParam(), and hitForFibre().
00338 { 00339 00340 G4StepPoint* preStepPoint = aStep->GetPreStepPoint(); 00341 const G4VTouchable* touch = preStepPoint->GetTouchable(); 00342 G4ThreeVector hitPoint = preStepPoint->GetPosition(); 00343 00344 int depth = (touch->GetReplicaNumber(0))%10 + 1; 00345 int lay = (touch->GetReplicaNumber(0)/10)%100 + 1; 00346 int det = (touch->GetReplicaNumber(1))/1000; 00347 00348 return setDetUnitId (det, hitPoint, depth, lay); 00349 }
void HCalSD::setNumberingScheme | ( | HcalNumberingScheme * | scheme | ) |
Definition at line 351 of file HCalSD.cc.
References numberingScheme.
Referenced by HCalSD(), HcalTB04Analysis::update(), SimG4HcalValidation::update(), and HcalTestAnalysis::update().
00351 { 00352 if (scheme != 0) { 00353 edm::LogInfo("HcalSim") << "HCalSD: updates numbering scheme for " 00354 << GetName() << "\n"; 00355 if (numberingScheme) delete numberingScheme; 00356 numberingScheme = scheme; 00357 } 00358 }
int HCalSD::setTrackID | ( | G4Step * | step | ) | [private] |
Definition at line 681 of file HCalSD.cc.
References TrackInformation::getIDonCaloSurface(), LogDebug, CaloSD::preStepPoint, CaloSD::previousID, CaloSD::resetForNewPrimary(), CaloSD::theTrack, and CaloHitID::trackID().
Referenced by getFromParam(), and hitForFibre().
00681 { 00682 00683 theTrack = aStep->GetTrack(); 00684 00685 double etrack = preStepPoint->GetKineticEnergy(); 00686 TrackInformation * trkInfo = (TrackInformation *)(theTrack->GetUserInformation()); 00687 int primaryID = trkInfo->getIDonCaloSurface(); 00688 if (primaryID == 0) { 00689 LogDebug("HcalSim") << "HCalSD: Problem with primaryID **** set by force " 00690 << "to TkID **** " << theTrack->GetTrackID(); 00691 primaryID = theTrack->GetTrackID(); 00692 } 00693 00694 if (primaryID != previousID.trackID()) 00695 resetForNewPrimary(preStepPoint->GetPosition(), etrack); 00696 00697 return primaryID; 00698 }
double HCalSD::betaThr [private] |
double HCalSD::birk1 [private] |
double HCalSD::birk2 [private] |
double HCalSD::birk3 [private] |
std::vector<G4LogicalVolume*> HCalSD::fibreLV [private] |
std::vector<G4String> HCalSD::fibreNames [private] |
std::vector<int> HCalSD::hfLevels [private] |
std::vector<G4LogicalVolume*> HCalSD::hfLV [private] |
std::vector<G4String> HCalSD::hfNames [private] |
HFShower* HCalSD::hfshower [private] |
std::vector<double> HCalSD::layer0wt [private] |
std::map<uint32_t,double> HCalSD::layerWeights [private] |
std::vector<G4Material*> HCalSD::materials [private] |
std::vector<G4String> HCalSD::matNames [private] |
G4int HCalSD::mumPDG [private] |
G4int HCalSD::mupPDG [private] |
HcalNumberingFromDDD* HCalSD::numberingFromDDD [private] |
Definition at line 65 of file HCalSD.h.
Referenced by getHitPMT(), HCalSD(), layerWeight(), setDetUnitId(), and ~HCalSD().
HcalNumberingScheme* HCalSD::numberingScheme [private] |
Definition at line 66 of file HCalSD.h.
Referenced by getHitPMT(), setDetUnitId(), setNumberingScheme(), and ~HCalSD().
std::vector<G4LogicalVolume*> HCalSD::pmtLV [private] |
std::vector<G4String> HCalSD::pmtNames [private] |
HFShowerLibrary* HCalSD::showerLibrary [private] |
HFShowerParam* HCalSD::showerParam [private] |
HFShowerPMT* HCalSD::showerPMT [private] |
Definition at line 70 of file HCalSD.h.
Referenced by getHitPMT(), HCalSD(), ProcessHits(), and ~HCalSD().
bool HCalSD::useBirk [private] |
bool HCalSD::useHF [private] |
bool HCalSD::useLayerWt [private] |
Definition at line 71 of file HCalSD.h.
Referenced by getEnergyDeposit(), HCalSD(), and readWeightFromFile().
bool HCalSD::useParam [private] |
bool HCalSD::usePMTHit [private] |
bool HCalSD::useShowerLibrary [private] |