CMS 3D CMS Logo

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