CMS 3D CMS Logo

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