CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC4_patch1/src/SimG4CMS/Calo/src/ECalSD.cc

Go to the documentation of this file.
00001 
00002 // File: ECalSD.cc
00003 // Description: Sensitive Detector class for electromagnetic calorimeters
00005 #include "SimG4CMS/Calo/interface/ECalSD.h"
00006 #include "SimG4Core/Notification/interface/TrackInformation.h"
00007 #include "Geometry/EcalCommonData/interface/EcalBarrelNumberingScheme.h"
00008 #include "Geometry/EcalCommonData/interface/EcalBaseNumber.h"
00009 #include "Geometry/EcalCommonData/interface/EcalEndcapNumberingScheme.h"
00010 #include "Geometry/EcalCommonData/interface/EcalPreshowerNumberingScheme.h"
00011 #include "Geometry/EcalCommonData/interface/ESTBNumberingScheme.h"
00012 #include "DataFormats/EcalDetId/interface/EBDetId.h"
00013 #include "DetectorDescription/Core/interface/DDFilter.h"
00014 #include "DetectorDescription/Core/interface/DDFilteredView.h"
00015 #include "DetectorDescription/Core/interface/DDSolid.h"
00016 #include "DetectorDescription/Core/interface/DDMaterial.h"
00017 #include "DetectorDescription/Core/interface/DDValue.h"
00018 
00019 #include "Geometry/EcalCommonData/interface/EcalBaseNumber.h"
00020 
00021 #include "G4LogicalVolumeStore.hh"
00022 #include "G4LogicalVolume.hh"
00023 #include "G4Step.hh"
00024 #include "G4Track.hh"
00025 #include "G4VProcess.hh"
00026 
00027 #include<algorithm>
00028 
00029 //#define DebugLog
00030 
00031 template <class T>
00032 bool any(const std::vector<T> & v, const T &what)
00033 {
00034   return std::find(v.begin(), v.end(), what) != v.end();
00035 }
00036 
00037 ECalSD::ECalSD(G4String name, const DDCompactView & cpv,
00038                SensitiveDetectorCatalog & clg, 
00039                edm::ParameterSet const & p, const SimTrackManager* manager) : 
00040   CaloSD(name, cpv, clg, p, manager, 
00041          p.getParameter<edm::ParameterSet>("ECalSD").getParameter<int>("TimeSliceUnit"),
00042          p.getParameter<edm::ParameterSet>("ECalSD").getParameter<bool>("IgnoreTrackID")), 
00043   numberingScheme(0){
00044 
00045   //   static SimpleConfigurable<bool>   on1(false,  "ECalSD:UseBirkLaw");
00046   //   static SimpleConfigurable<double> bk1(0.00463,"ECalSD:BirkC1");
00047   //   static SimpleConfigurable<double> bk2(-0.03,  "ECalSD:BirkC2");
00048   //   static SimpleConfigurable<double> bk3(1.0,    "ECalSD:BirkC3");
00049   // Values from NIM A484 (2002) 239-244: as implemented in Geant3
00050   //   useBirk          = on1.value();
00051   //   birk1            = bk1.value()*(g/(MeV*cm2));
00052   //   birk2            = bk2.value()*(g/(MeV*cm2))*(g/(MeV*cm2));
00053 
00054   edm::ParameterSet m_EC = p.getParameter<edm::ParameterSet>("ECalSD");
00055   useBirk      = m_EC.getParameter<bool>("UseBirkLaw");
00056   useBirkL3    = m_EC.getParameter<bool>("BirkL3Parametrization");
00057   birk1        = m_EC.getParameter<double>("BirkC1")*(g/(MeV*cm2));
00058   birk2        = m_EC.getParameter<double>("BirkC2");
00059   birk3        = m_EC.getParameter<double>("BirkC3");
00060   birkSlope    = m_EC.getParameter<double>("BirkSlope");
00061   birkCut      = m_EC.getParameter<double>("BirkCut");
00062   slopeLY      = m_EC.getParameter<double>("SlopeLightYield");
00063   storeTrack   = m_EC.getParameter<bool>("StoreSecondary");
00064   crystalMat   = m_EC.getUntrackedParameter<std::string>("XtalMat","E_PbWO4");
00065   bool isItTB  = m_EC.getUntrackedParameter<bool>("TestBeam", false);
00066   bool nullNS  = m_EC.getUntrackedParameter<bool>("NullNumbering", false);
00067   storeRL      = m_EC.getUntrackedParameter<bool>("StoreRadLength", false);
00068   
00069   ageingWithSlopeLY   = m_EC.getUntrackedParameter<bool>("AgeingWithSlopeLY", false);
00070   if(ageingWithSlopeLY) ageing.setLumies(p.getParameter<edm::ParameterSet>("ECalSD").getParameter<double>("DelivLuminosity"),
00071                                          p.getParameter<edm::ParameterSet>("ECalSD").getParameter<double>("InstLuminosity"));
00072 
00073   //Material list for HB/HE/HO sensitive detectors
00074   std::string attribute = "ReadOutName";
00075   DDSpecificsFilter filter;
00076   DDValue           ddv(attribute,name,0);
00077   filter.setCriteria(ddv,DDSpecificsFilter::equals);
00078   DDFilteredView fv(cpv);
00079   fv.addFilter(filter);
00080   fv.firstChild();
00081   DDsvalues_type sv(fv.mergedSpecifics());
00082   // Use of Weight
00083   useWeight= true;
00084   std::vector<double>      tempD = getDDDArray("EnergyWeight",sv);
00085   if (tempD.size() > 0) { if (tempD[0] < 0.1) useWeight = false; }
00086   std::vector<std::string> tempS = getStringArray("Depth1Name",sv);
00087   if (tempS.size() > 0) depth1Name = tempS[0];
00088   else                  depth1Name = " ";
00089   tempS = getStringArray("Depth2Name",sv);
00090   if (tempS.size() > 0) depth2Name = tempS[0];
00091   else                  depth2Name = " ";
00092 
00093   EcalNumberingScheme* scheme=0;
00094   if (nullNS)                    scheme = 0;
00095   else if (name == "EcalHitsEB") scheme = dynamic_cast<EcalNumberingScheme*>(new EcalBarrelNumberingScheme());
00096   else if (name == "EcalHitsEE") scheme = dynamic_cast<EcalNumberingScheme*>(new EcalEndcapNumberingScheme());
00097   else if (name == "EcalHitsES") {
00098     if (isItTB) scheme = dynamic_cast<EcalNumberingScheme*>(new ESTBNumberingScheme());
00099     else        scheme = dynamic_cast<EcalNumberingScheme*>(new EcalPreshowerNumberingScheme());
00100     useWeight = false;
00101   } else {edm::LogWarning("EcalSim") << "ECalSD: ReadoutName not supported\n";}
00102 
00103   if (scheme)  setNumberingScheme(scheme);
00104 #ifdef DebugLog
00105   LogDebug("EcalSim") << "Constructing a ECalSD  with name " << GetName();
00106 #endif
00107   if (useWeight) {
00108     edm::LogInfo("EcalSim")  << "ECalSD:: Use of Birks law is set to      " 
00109                              << useBirk << "        with three constants kB = "
00110                              << birk1 << ", C1 = " << birk2 << ", C2 = " 
00111                              << birk3 <<"\n         Use of L3 parametrization "
00112                              << useBirkL3 << " with slope " << birkSlope
00113                              << " and cut off " << birkCut << "\n"
00114                              << "         Slope for Light yield is set to "
00115                              << slopeLY;
00116   } else {
00117     edm::LogInfo("EcalSim")  << "ECalSD:: energy deposit is not corrected "
00118                              << " by Birk or light yield curve";
00119   }
00120 
00121   edm::LogInfo("EcalSim") << "ECalSD:: Suppression Flag " << suppressHeavy
00122                           << " protons below " << kmaxProton << " MeV,"
00123                           << " neutrons below " << kmaxNeutron << " MeV and"
00124                           << " ions below " << kmaxIon << " MeV\n"
00125                           << "         Depth1 Name = " << depth1Name
00126                           << " and Depth2 Name = " << depth2Name;
00127 
00128   if (useWeight) initMap(name,cpv);
00129 
00130 }
00131 
00132 ECalSD::~ECalSD() {
00133   if (numberingScheme) delete numberingScheme;
00134 }
00135 
00136 double ECalSD::getEnergyDeposit(G4Step * aStep) {
00137   
00138   if (aStep == NULL) {
00139     return 0;
00140   } else {
00141     preStepPoint        = aStep->GetPreStepPoint();
00142     G4Track* theTrack   = aStep->GetTrack();
00143     double wt2  = theTrack->GetWeight();
00144     G4String nameVolume = preStepPoint->GetPhysicalVolume()->GetName();
00145 
00146     // take into account light collection curve for crystals
00147     double weight = 1.;
00148     if (suppressHeavy) {
00149       TrackInformation * trkInfo = (TrackInformation *)(theTrack->GetUserInformation());
00150       if (trkInfo) {
00151         int pdg = theTrack->GetDefinition()->GetPDGEncoding();
00152         if (!(trkInfo->isPrimary())) { // Only secondary particles
00153           double ke = theTrack->GetKineticEnergy()/MeV;
00154           if (((pdg/1000000000 == 1 && ((pdg/10000)%100) > 0 &&
00155                 ((pdg/10)%100) > 0)) && (ke<kmaxIon)) weight = 0;
00156           if ((pdg == 2212) && (ke < kmaxProton))     weight = 0;
00157           if ((pdg == 2112) && (ke < kmaxNeutron))    weight = 0;
00158 #ifdef DebugLog
00159           if (weight == 0) 
00160             LogDebug("EcalSim") << "Ignore Track " << theTrack->GetTrackID()
00161                                 << " Type " << theTrack->GetDefinition()->GetParticleName()
00162                                 << " Kinetic Energy " << ke << " MeV";
00163 #endif
00164         }
00165       }
00166     }
00167     G4LogicalVolume* lv   = aStep->GetPreStepPoint()->GetTouchable()->GetVolume(0)->GetLogicalVolume();
00168     if (useWeight && !any(noWeight,lv)) {
00169       weight *= curve_LY(aStep);
00170       if (useBirk) {
00171         if (useBirkL3) weight *= getBirkL3(aStep);
00172         else           weight *= getAttenuation(aStep, birk1, birk2, birk3);
00173       }
00174     }
00175     double wt1  = getResponseWt(theTrack);
00176     double edep = aStep->GetTotalEnergyDeposit()*weight*wt1;
00177     /*
00178     if(wt2 != 1.0) { 
00179       std::cout << "ECalSD:: " << nameVolume
00180                 <<" LightWeight= " <<weight << " wt1= " <<wt1
00181                 << "  wt2= " << wt2 << "  "
00182                 << " Weighted Energy Deposit " << edep/MeV << " MeV" 
00183                 << std::endl;
00184       std::cout << theTrack->GetDefinition()->GetParticleName()
00185                 << " " << theTrack->GetKineticEnergy()
00186                 << " Id=" << theTrack->GetTrackID()
00187                 << " IdP=" << theTrack->GetParentID();
00188       const G4VProcess* pr = theTrack->GetCreatorProcess();
00189       if(pr) std::cout << " from  " << pr->GetProcessName();
00190       std::cout << std::endl;
00191     }
00192     */
00193     if(wt2 > 0.0) { edep *= wt2; }
00194 #ifdef DebugLog
00195     LogDebug("EcalSim") << "ECalSD:: " << nameVolume
00196                         <<" Light Collection Efficiency " <<weight << ":" <<wt1
00197                         << " Weighted Energy Deposit " << edep/MeV << " MeV";
00198 #endif
00199     return edep;
00200   } 
00201 }
00202 
00203 int ECalSD::getTrackID(G4Track* aTrack) {
00204 
00205   int  primaryID(0);
00206   bool flag(false);
00207   if (storeTrack) {
00208     G4LogicalVolume* lv  = preStepPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
00209     if (any(useDepth1,lv)) {
00210       flag = true;
00211     } else if (any(useDepth2,lv)) {
00212       flag = true;
00213     }
00214   }
00215   if (flag) {
00216     forceSave = true;
00217     primaryID = aTrack->GetTrackID();
00218   } else {
00219     primaryID = CaloSD::getTrackID(aTrack);
00220   }
00221   return primaryID;
00222 }
00223 
00224 uint16_t ECalSD::getDepth(G4Step * aStep) {
00225   G4LogicalVolume* lv   = aStep->GetPreStepPoint()->GetTouchable()->GetVolume(0)->GetLogicalVolume();
00226   uint16_t ret = 0;
00227   if (any(useDepth1,lv))      ret = 1;
00228   else if (any(useDepth2,lv)) ret = 2;
00229   else if (storeRL) ret = getRadiationLength(aStep);
00230 #ifdef DebugLog
00231   LogDebug("EcalSim") << "Volume " << lv->GetName() << " Depth " << ret;
00232 #endif
00233   return ret;
00234 }
00235 
00236 uint16_t ECalSD::getRadiationLength(G4Step * aStep) {
00237   
00238   uint16_t thisX0 = 0;
00239   if (aStep != NULL) {
00240     G4StepPoint* hitPoint = aStep->GetPreStepPoint();
00241     G4LogicalVolume* lv   = hitPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
00242     
00243     if (useWeight) {
00244       G4ThreeVector  localPoint = setToLocal(hitPoint->GetPosition(),
00245                                              hitPoint->GetTouchable());
00246       double crlength = crystalLength(lv);
00247       double radl     = hitPoint->GetMaterial()->GetRadlen();
00248       double detz     = (float)(0.5*crlength + localPoint.z());
00249       thisX0 = (uint16_t)floor(detz/radl);   
00250     } 
00251   }
00252   return thisX0;
00253 }
00254 
00255 uint32_t ECalSD::setDetUnitId(G4Step * aStep) { 
00256   if (numberingScheme == 0) {
00257     return EBDetId(1,1)();
00258   } else {
00259     getBaseNumber(aStep);
00260     return numberingScheme->getUnitID(theBaseNumber);
00261   }
00262 }
00263 
00264 void ECalSD::setNumberingScheme(EcalNumberingScheme* scheme) {
00265   if (scheme != 0) {
00266     edm::LogInfo("EcalSim") << "EcalSD: updates numbering scheme for " 
00267                             << GetName() << "\n";
00268     if (numberingScheme) delete numberingScheme;
00269     numberingScheme = scheme;
00270   }
00271 }
00272 
00273 
00274 void ECalSD::initMap(G4String sd, const DDCompactView & cpv) {
00275 
00276   G4String attribute = "ReadOutName";
00277   DDSpecificsFilter filter;
00278   DDValue           ddv(attribute,sd,0);
00279   filter.setCriteria(ddv,DDSpecificsFilter::equals);
00280   DDFilteredView fv(cpv);
00281   fv.addFilter(filter);
00282   fv.firstChild();
00283 
00284   std::vector<G4LogicalVolume*> lvused;
00285   const G4LogicalVolumeStore *  lvs = G4LogicalVolumeStore::GetInstance();
00286   std::vector<G4LogicalVolume *>::const_iterator lvcite;
00287   std::map<std::string, G4LogicalVolume *> nameMap;
00288   for (auto lvi = lvs->begin(), lve = lvs->end(); lvi != lve; ++lvi)
00289     nameMap.insert(std::make_pair((*lvi)->GetName(), *lvi));
00290 
00291   bool dodet=true;
00292   while (dodet) {
00293     const std::string &matname = fv.logicalPart().material().name().name();
00294     const std::string &lvname = fv.logicalPart().name().name();
00295     G4LogicalVolume* lv = nameMap[lvname];
00296     if (depth1Name != " ") {
00297       if (strncmp(lvname.c_str(), depth1Name.c_str(), 4) == 0) {
00298         if (!any(useDepth1, lv)) {
00299           useDepth1.push_back(lv);
00300 #ifdef DebugLog
00301           LogDebug("EcalSim") << "ECalSD::initMap Logical Volume " << lvname
00302                               <<" in Depth 1 volume list";
00303 #endif
00304         }
00305         G4LogicalVolume* lvr = nameMap[lvname + "_refl"];
00306         if (lvr != 0 && !any(useDepth1, lvr)) {
00307           useDepth1.push_back(lvr);
00308 #ifdef DebugLog
00309           LogDebug("EcalSim") << "ECalSD::initMap Logical Volume " << lvname << "_refl"
00310                               <<" in Depth 1 volume list";
00311 #endif
00312         }
00313       }
00314     }
00315     if (depth2Name != " ") {
00316       if (strncmp(lvname.c_str(), depth2Name.c_str(), 4) == 0) {
00317         if (!any(useDepth2, lv)) {
00318           useDepth2.push_back(lv);
00319 #ifdef DebugLog
00320           LogDebug("EcalSim") << "ECalSD::initMap Logical Volume " << lvname
00321                               <<" in Depth 2 volume list";
00322 #endif
00323         }
00324         G4LogicalVolume* lvr = nameMap[lvname + "_refl"];
00325         if (lvr != 0 && !any(useDepth2,lvr)) {
00326           useDepth2.push_back(lvr);
00327 #ifdef DebugLog
00328           LogDebug("EcalSim") << "ECalSD::initMap Logical Volume " << lvname << "_refl"
00329                               <<" in Depth 2 volume list";
00330 #endif
00331         }
00332       }
00333     }
00334     if (lv != 0) {
00335       if (crystalMat.size() == matname.size() && !strcmp(crystalMat.c_str(), matname.c_str())) {
00336         if (!any(lvused,lv)) {
00337           lvused.push_back(lv);
00338           const DDSolid & sol  = fv.logicalPart().solid();
00339           const std::vector<double> & paras = sol.parameters();
00340 #ifdef DebugLog
00341           LogDebug("EcalSim") << "ECalSD::initMap (for " << sd << "): Solid " 
00342                               << lvname << " Shape " << sol.shape() 
00343                               << " Parameter 0 = "<< paras[0] 
00344                               << " Logical Volume " << lv;
00345 #endif
00346           if (sol.shape() == ddtrap) {
00347             double dz = 2*paras[0];
00348             xtalLMap.insert(std::pair<G4LogicalVolume*,double>(lv,dz));
00349             lv = nameMap[lvname + "_refl"];
00350             if (lv != 0)
00351               xtalLMap.insert(std::pair<G4LogicalVolume*,double>(lv,dz));
00352           }
00353         }
00354       } else {
00355         if (!any(noWeight,lv)) {
00356           noWeight.push_back(lv);
00357 #ifdef DebugLog
00358           LogDebug("EcalSim") << "ECalSD::initMap Logical Volume " << lvname
00359                               << " Material " << matname <<" in noWeight list";
00360 #endif
00361         }
00362         lv = nameMap[lvname];
00363         if (lv != 0 && !any(noWeight,lv)) {
00364           noWeight.push_back(lv);
00365 #ifdef DebugLog
00366           LogDebug("EcalSim") << "ECalSD::initMap Logical Volume " << lvname
00367                               << " Material " << matname <<" in noWeight list";
00368 #endif
00369         }
00370       }
00371     }
00372     dodet = fv.next();
00373   }
00374 #ifdef DebugLog
00375   LogDebug("EcalSim") << "ECalSD: Length Table for " << attribute << " = " 
00376                       << sd << ":";   
00377   std::map<G4LogicalVolume*,double>::const_iterator ite = xtalLMap.begin();
00378   int i=0;
00379   for (; ite != xtalLMap.end(); ite++, i++) {
00380     G4String name = "Unknown";
00381     if (ite->first != 0) name = (ite->first)->GetName();
00382     LogDebug("EcalSim") << " " << i << " " << ite->first << " " << name 
00383                         << " L = " << ite->second;
00384   }
00385 #endif
00386 }
00387 
00388 double ECalSD::curve_LY(G4Step* aStep) {
00389 
00390   G4StepPoint*     stepPoint = aStep->GetPreStepPoint();
00391   G4LogicalVolume* lv        = stepPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
00392 
00393   double weight = 1.;
00394   G4ThreeVector  localPoint = setToLocal(stepPoint->GetPosition(),
00395                                          stepPoint->GetTouchable());
00396 
00397   double crlength = crystalLength(lv);
00398 
00399   if(ageingWithSlopeLY){
00400     //position along the crystal in mm from 0 to 230 (in EB)
00401     double depth = 0.5 * crlength + localPoint.z();
00402 
00403     if (depth >= -0.1 || depth <= crlength+0.1)
00404       weight = ageing.calcLightCollectionEfficiencyWeighted(currentID.unitID(), depth/crlength);
00405   }
00406   else{
00407     double dapd = 0.5 * crlength - localPoint.z();
00408     if (dapd >= -0.1 || dapd <= crlength+0.1) {
00409       if (dapd <= 100.)
00410         weight = 1.0 + slopeLY - dapd * 0.01 * slopeLY;
00411     } else {
00412       edm::LogWarning("EcalSim") << "ECalSD: light coll curve : wrong distance "
00413                                  << "to APD " << dapd << " crlength = " 
00414                                  << crlength <<" crystal name = " <<lv->GetName()
00415                                  << " z of localPoint = " << localPoint.z() 
00416                                  << " take weight = " << weight;
00417     }
00418   }
00419 #ifdef DebugLog
00420   LogDebug("EcalSim") << "ECalSD, light coll curve : " << dapd 
00421                       << " crlength = " << crlength
00422                       << " crystal name = " << lv->GetName()
00423                       << " z of localPoint = " << localPoint.z() 
00424                       << " take weight = " << weight;
00425 #endif
00426   return weight;
00427 }
00428 
00429 double ECalSD::crystalLength(G4LogicalVolume* lv) {
00430 
00431   double length= 230.;
00432   std::map<G4LogicalVolume*,double>::const_iterator ite = xtalLMap.find(lv);
00433   if (ite != xtalLMap.end()) length = ite->second;
00434   return length;
00435 }
00436 
00437 void ECalSD::getBaseNumber(const G4Step* aStep) {
00438 
00439   theBaseNumber.reset();
00440   const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
00441   int theSize = touch->GetHistoryDepth()+1;
00442   if ( theBaseNumber.getCapacity() < theSize ) theBaseNumber.setSize(theSize);
00443   //Get name and copy numbers
00444   if ( theSize > 1 ) {
00445     for (int ii = 0; ii < theSize ; ii++) {
00446       theBaseNumber.addLevel(touch->GetVolume(ii)->GetName(),touch->GetReplicaNumber(ii));
00447 #ifdef DebugLog
00448       LogDebug("EcalSim") << "ECalSD::getBaseNumber(): Adding level " << ii
00449                           << ": " << touch->GetVolume(ii)->GetName() << "["
00450                           << touch->GetReplicaNumber(ii) << "]";
00451 #endif
00452     }
00453   }
00454 
00455 }
00456 
00457 double ECalSD::getBirkL3(G4Step* aStep) {
00458 
00459   double weight = 1.;
00460   double charge = aStep->GetPreStepPoint()->GetCharge();
00461 
00462   if (charge != 0. && aStep->GetStepLength() > 0) {
00463     G4Material* mat = aStep->GetPreStepPoint()->GetMaterial();
00464     double density = mat->GetDensity();
00465     double dedx    = aStep->GetTotalEnergyDeposit()/aStep->GetStepLength();
00466     double rkb     = birk1/density;
00467     if (dedx > 0) {
00468       weight         = 1. - birkSlope*log(rkb*dedx);
00469       if (weight < birkCut) weight = birkCut;
00470       else if (weight > 1.) weight = 1.;
00471     }
00472 #ifdef DebugLog
00473     LogDebug("EcalSim") << "ECalSD::getBirkL3 in " << mat->GetName()
00474                         << " Charge " << charge << " dE/dx " << dedx
00475                         << " Birk Const " << rkb << " Weight = " << weight 
00476                         << " dE " << aStep->GetTotalEnergyDeposit();
00477 #endif
00478   }
00479   return weight;
00480 
00481 }
00482 
00483 std::vector<double> ECalSD::getDDDArray(const std::string & str,
00484                                         const DDsvalues_type & sv) {
00485 
00486 #ifdef DebugLog
00487   LogDebug("EcalSim") << "ECalSD:getDDDArray called for " << str;
00488 #endif
00489   DDValue value(str);
00490   if (DDfetch(&sv,value)) {
00491 #ifdef DebugLog
00492     LogDebug("EcalSim") << value;
00493 #endif
00494     const std::vector<double> & fvec = value.doubles();
00495     return fvec;
00496   } else {
00497     std::vector<double> fvec;
00498     return fvec;
00499   }
00500 }
00501 
00502 std::vector<std::string> ECalSD::getStringArray(const std::string & str,
00503                                                 const DDsvalues_type & sv) {
00504 
00505 #ifdef DebugLog
00506   LogDebug("EcalSim") << "ECalSD:getStringArray called for " << str;
00507 #endif
00508   DDValue value(str);
00509   if (DDfetch(&sv,value)) {
00510 #ifdef DebugLog
00511     LogDebug("EcalSim") << value;
00512 #endif
00513     const std::vector<std::string> & fvec = value.strings();
00514     return fvec;
00515   } else {
00516     std::vector<std::string> fvec;
00517     return fvec;
00518   }
00519 }