CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_1/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     G4Track* theTrack   = aStep->GetTrack();
00139     double wt2  = theTrack->GetWeight();
00140     G4String nameVolume = preStepPoint->GetPhysicalVolume()->GetName();
00141 
00142     // take into account light collection curve for crystals
00143     double weight = 1.;
00144     if (suppressHeavy) {
00145       TrackInformation * trkInfo = (TrackInformation *)(theTrack->GetUserInformation());
00146       if (trkInfo) {
00147         int pdg = theTrack->GetDefinition()->GetPDGEncoding();
00148         if (!(trkInfo->isPrimary())) { // Only secondary particles
00149           double ke = theTrack->GetKineticEnergy()/MeV;
00150           if (((pdg/1000000000 == 1 && ((pdg/10000)%100) > 0 &&
00151                 ((pdg/10)%100) > 0)) && (ke<kmaxIon)) weight = 0;
00152           if ((pdg == 2212) && (ke < kmaxProton))     weight = 0;
00153           if ((pdg == 2112) && (ke < kmaxNeutron))    weight = 0;
00154 #ifdef DebugLog
00155           if (weight == 0) 
00156             LogDebug("EcalSim") << "Ignore Track " << theTrack->GetTrackID()
00157                                 << " Type " << theTrack->GetDefinition()->GetParticleName()
00158                                 << " Kinetic Energy " << ke << " MeV";
00159 #endif
00160         }
00161       }
00162     }
00163     G4LogicalVolume* lv   = aStep->GetPreStepPoint()->GetTouchable()->GetVolume(0)->GetLogicalVolume();
00164     if (useWeight && !any(noWeight,lv)) {
00165       weight *= curve_LY(aStep);
00166       if (useBirk) {
00167         if (useBirkL3) weight *= getBirkL3(aStep);
00168         else           weight *= getAttenuation(aStep, birk1, birk2, birk3);
00169       }
00170     }
00171     double wt1  = getResponseWt(theTrack);
00172     double edep = aStep->GetTotalEnergyDeposit()*weight*wt1;
00173     /*
00174     if(wt2 != 1.0) { 
00175       std::cout << "ECalSD:: " << nameVolume
00176                 <<" LightWeight= " <<weight << " wt1= " <<wt1
00177                 << "  wt2= " << wt2 << "  "
00178                 << " Weighted Energy Deposit " << edep/MeV << " MeV" 
00179                 << std::endl;
00180       std::cout << theTrack->GetDefinition()->GetParticleName()
00181                 << " " << theTrack->GetKineticEnergy()
00182                 << " Id=" << theTrack->GetTrackID()
00183                 << " IdP=" << theTrack->GetParentID();
00184       const G4VProcess* pr = theTrack->GetCreatorProcess();
00185       if(pr) std::cout << " from  " << pr->GetProcessName();
00186       std::cout << std::endl;
00187     }
00188     */
00189     if(wt2 > 0.0) { edep *= wt2; }
00190 #ifdef DebugLog
00191     LogDebug("EcalSim") << "ECalSD:: " << nameVolume
00192                         <<" Light Collection Efficiency " <<weight << ":" <<wt1
00193                         << " Weighted Energy Deposit " << edep/MeV << " MeV";
00194 #endif
00195     return edep;
00196   } 
00197 }
00198 
00199 int ECalSD::getTrackID(G4Track* aTrack) {
00200 
00201   int  primaryID(0);
00202   bool flag(false);
00203   if (storeTrack) {
00204     G4LogicalVolume* lv  = preStepPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
00205     if (any(useDepth1,lv)) {
00206       flag = true;
00207     } else if (any(useDepth2,lv)) {
00208       flag = true;
00209     }
00210   }
00211   if (flag) {
00212     forceSave = true;
00213     primaryID = aTrack->GetTrackID();
00214   } else {
00215     primaryID = CaloSD::getTrackID(aTrack);
00216   }
00217   return primaryID;
00218 }
00219 
00220 uint16_t ECalSD::getDepth(G4Step * aStep) {
00221   G4LogicalVolume* lv   = aStep->GetPreStepPoint()->GetTouchable()->GetVolume(0)->GetLogicalVolume();
00222   uint16_t ret = 0;
00223   if (any(useDepth1,lv))      ret = 1;
00224   else if (any(useDepth2,lv)) ret = 2;
00225   else if (storeRL) ret = getRadiationLength(aStep);
00226 #ifdef DebugLog
00227   LogDebug("EcalSim") << "Volume " << lv->GetName() << " Depth " << ret;
00228 #endif
00229   return ret;
00230 }
00231 
00232 uint16_t ECalSD::getRadiationLength(G4Step * aStep) {
00233   
00234   uint16_t thisX0 = 0;
00235   if (aStep != NULL) {
00236     G4StepPoint* hitPoint = aStep->GetPreStepPoint();
00237     G4LogicalVolume* lv   = hitPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
00238     
00239     if (useWeight) {
00240       G4ThreeVector  localPoint = setToLocal(hitPoint->GetPosition(),
00241                                              hitPoint->GetTouchable());
00242       double crlength = crystalLength(lv);
00243       double radl     = hitPoint->GetMaterial()->GetRadlen();
00244       double detz     = (float)(0.5*crlength + localPoint.z());
00245       thisX0 = (uint16_t)floor(detz/radl);   
00246     } 
00247   }
00248   return thisX0;
00249 }
00250 
00251 uint32_t ECalSD::setDetUnitId(G4Step * aStep) { 
00252   if (numberingScheme == 0) {
00253     return EBDetId(1,1)();
00254   } else {
00255     getBaseNumber(aStep);
00256     return numberingScheme->getUnitID(theBaseNumber);
00257   }
00258 }
00259 
00260 void ECalSD::setNumberingScheme(EcalNumberingScheme* scheme) {
00261   if (scheme != 0) {
00262     edm::LogInfo("EcalSim") << "EcalSD: updates numbering scheme for " 
00263                             << GetName() << "\n";
00264     if (numberingScheme) delete numberingScheme;
00265     numberingScheme = scheme;
00266   }
00267 }
00268 
00269 
00270 void ECalSD::initMap(G4String sd, const DDCompactView & cpv) {
00271 
00272   G4String attribute = "ReadOutName";
00273   DDSpecificsFilter filter;
00274   DDValue           ddv(attribute,sd,0);
00275   filter.setCriteria(ddv,DDSpecificsFilter::equals);
00276   DDFilteredView fv(cpv);
00277   fv.addFilter(filter);
00278   fv.firstChild();
00279 
00280   std::vector<G4LogicalVolume*> lvused;
00281   const G4LogicalVolumeStore *  lvs = G4LogicalVolumeStore::GetInstance();
00282   std::vector<G4LogicalVolume *>::const_iterator lvcite;
00283   std::map<std::string, G4LogicalVolume *> nameMap;
00284   for (auto lvi = lvs->begin(), lve = lvs->end(); lvi != lve; ++lvi)
00285     nameMap.insert(std::make_pair((*lvi)->GetName(), *lvi));
00286 
00287   bool dodet=true;
00288   while (dodet) {
00289     const std::string &matname = fv.logicalPart().material().name().name();
00290     const std::string &lvname = fv.logicalPart().name().name();
00291     G4LogicalVolume* lv = nameMap[lvname];
00292     if (depth1Name != " ") {
00293       if (strncmp(lvname.c_str(), depth1Name.c_str(), 4) == 0) {
00294         if (!any(useDepth1, lv)) {
00295           useDepth1.push_back(lv);
00296 #ifdef DebugLog
00297           LogDebug("EcalSim") << "ECalSD::initMap Logical Volume " << lvname
00298                               <<" in Depth 1 volume list";
00299 #endif
00300         }
00301         G4LogicalVolume* lvr = nameMap[lvname + "_refl"];
00302         if (lvr != 0 && !any(useDepth1, lvr)) {
00303           useDepth1.push_back(lvr);
00304 #ifdef DebugLog
00305           LogDebug("EcalSim") << "ECalSD::initMap Logical Volume " << lvname << "_refl"
00306                               <<" in Depth 1 volume list";
00307 #endif
00308         }
00309       }
00310     }
00311     if (depth2Name != " ") {
00312       if (strncmp(lvname.c_str(), depth2Name.c_str(), 4) == 0) {
00313         if (!any(useDepth2, lv)) {
00314           useDepth2.push_back(lv);
00315 #ifdef DebugLog
00316           LogDebug("EcalSim") << "ECalSD::initMap Logical Volume " << lvname
00317                               <<" in Depth 2 volume list";
00318 #endif
00319         }
00320         G4LogicalVolume* lvr = nameMap[lvname + "_refl"];
00321         if (lvr != 0 && !any(useDepth2,lvr)) {
00322           useDepth2.push_back(lvr);
00323 #ifdef DebugLog
00324           LogDebug("EcalSim") << "ECalSD::initMap Logical Volume " << lvname << "_refl"
00325                               <<" in Depth 2 volume list";
00326 #endif
00327         }
00328       }
00329     }
00330     if (lv != 0) {
00331       if (crystalMat.size() == matname.size() && !strcmp(crystalMat.c_str(), matname.c_str())) {
00332         if (!any(lvused,lv)) {
00333           lvused.push_back(lv);
00334           const DDSolid & sol  = fv.logicalPart().solid();
00335           const std::vector<double> & paras = sol.parameters();
00336 #ifdef DebugLog
00337           LogDebug("EcalSim") << "ECalSD::initMap (for " << sd << "): Solid " 
00338                               << lvname << " Shape " << sol.shape() 
00339                               << " Parameter 0 = "<< paras[0] 
00340                               << " Logical Volume " << lv;
00341 #endif
00342           if (sol.shape() == ddtrap) {
00343             double dz = 2*paras[0];
00344             xtalLMap.insert(std::pair<G4LogicalVolume*,double>(lv,dz));
00345             lv = nameMap[lvname + "_refl"];
00346             if (lv != 0)
00347               xtalLMap.insert(std::pair<G4LogicalVolume*,double>(lv,dz));
00348           }
00349         }
00350       } else {
00351         if (!any(noWeight,lv)) {
00352           noWeight.push_back(lv);
00353 #ifdef DebugLog
00354           LogDebug("EcalSim") << "ECalSD::initMap Logical Volume " << lvname
00355                               << " Material " << matname <<" in noWeight list";
00356 #endif
00357         }
00358         lv = nameMap[lvname];
00359         if (lv != 0 && !any(noWeight,lv)) {
00360           noWeight.push_back(lv);
00361 #ifdef DebugLog
00362           LogDebug("EcalSim") << "ECalSD::initMap Logical Volume " << lvname
00363                               << " Material " << matname <<" in noWeight list";
00364 #endif
00365         }
00366       }
00367     }
00368     dodet = fv.next();
00369   }
00370 #ifdef DebugLog
00371   LogDebug("EcalSim") << "ECalSD: Length Table for " << attribute << " = " 
00372                       << sd << ":";   
00373   std::map<G4LogicalVolume*,double>::const_iterator ite = xtalLMap.begin();
00374   int i=0;
00375   for (; ite != xtalLMap.end(); ite++, i++) {
00376     G4String name = "Unknown";
00377     if (ite->first != 0) name = (ite->first)->GetName();
00378     LogDebug("EcalSim") << " " << i << " " << ite->first << " " << name 
00379                         << " L = " << ite->second;
00380   }
00381 #endif
00382 }
00383 
00384 double ECalSD::curve_LY(G4Step* aStep) {
00385 
00386   G4StepPoint*     stepPoint = aStep->GetPreStepPoint();
00387   G4LogicalVolume* lv        = stepPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
00388 
00389   double weight = 1.;
00390   G4ThreeVector  localPoint = setToLocal(stepPoint->GetPosition(),
00391                                          stepPoint->GetTouchable());
00392   double crlength = crystalLength(lv);
00393   double dapd = 0.5 * crlength - localPoint.z();
00394   if (dapd >= -0.1 || dapd <= crlength+0.1) {
00395     if (dapd <= 100.)
00396       weight = 1.0 + slopeLY - dapd * 0.01 * slopeLY;
00397   } else {
00398     edm::LogWarning("EcalSim") << "ECalSD: light coll curve : wrong distance "
00399                                << "to APD " << dapd << " crlength = " 
00400                                << crlength <<" crystal name = " <<lv->GetName()
00401                                << " z of localPoint = " << localPoint.z() 
00402                                << " take weight = " << weight;
00403   }
00404 #ifdef DebugLog
00405   LogDebug("EcalSim") << "ECalSD, light coll curve : " << dapd 
00406                       << " crlength = " << crlength
00407                       << " crystal name = " << lv->GetName()
00408                       << " z of localPoint = " << localPoint.z() 
00409                       << " take weight = " << weight;
00410 #endif
00411   return weight;
00412 }
00413 
00414 double ECalSD::crystalLength(G4LogicalVolume* lv) {
00415 
00416   double length= 230.;
00417   std::map<G4LogicalVolume*,double>::const_iterator ite = xtalLMap.find(lv);
00418   if (ite != xtalLMap.end()) length = ite->second;
00419   return length;
00420 }
00421 
00422 void ECalSD::getBaseNumber(const G4Step* aStep) {
00423 
00424   theBaseNumber.reset();
00425   const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
00426   int theSize = touch->GetHistoryDepth()+1;
00427   if ( theBaseNumber.getCapacity() < theSize ) theBaseNumber.setSize(theSize);
00428   //Get name and copy numbers
00429   if ( theSize > 1 ) {
00430     for (int ii = 0; ii < theSize ; ii++) {
00431       theBaseNumber.addLevel(touch->GetVolume(ii)->GetName(),touch->GetReplicaNumber(ii));
00432 #ifdef DebugLog
00433       LogDebug("EcalSim") << "ECalSD::getBaseNumber(): Adding level " << ii
00434                           << ": " << touch->GetVolume(ii)->GetName() << "["
00435                           << touch->GetReplicaNumber(ii) << "]";
00436 #endif
00437     }
00438   }
00439 
00440 }
00441 
00442 double ECalSD::getBirkL3(G4Step* aStep) {
00443 
00444   double weight = 1.;
00445   double charge = aStep->GetPreStepPoint()->GetCharge();
00446 
00447   if (charge != 0. && aStep->GetStepLength() > 0) {
00448     G4Material* mat = aStep->GetPreStepPoint()->GetMaterial();
00449     double density = mat->GetDensity();
00450     double dedx    = aStep->GetTotalEnergyDeposit()/aStep->GetStepLength();
00451     double rkb     = birk1/density;
00452     if (dedx > 0) {
00453       weight         = 1. - birkSlope*log(rkb*dedx);
00454       if (weight < birkCut) weight = birkCut;
00455       else if (weight > 1.) weight = 1.;
00456     }
00457 #ifdef DebugLog
00458     LogDebug("EcalSim") << "ECalSD::getBirkL3 in " << mat->GetName()
00459                         << " Charge " << charge << " dE/dx " << dedx
00460                         << " Birk Const " << rkb << " Weight = " << weight 
00461                         << " dE " << aStep->GetTotalEnergyDeposit();
00462 #endif
00463   }
00464   return weight;
00465 
00466 }
00467 
00468 std::vector<double> ECalSD::getDDDArray(const std::string & str,
00469                                         const DDsvalues_type & sv) {
00470 
00471 #ifdef DebugLog
00472   LogDebug("EcalSim") << "ECalSD:getDDDArray called for " << str;
00473 #endif
00474   DDValue value(str);
00475   if (DDfetch(&sv,value)) {
00476 #ifdef DebugLog
00477     LogDebug("EcalSim") << value;
00478 #endif
00479     const std::vector<double> & fvec = value.doubles();
00480     return fvec;
00481   } else {
00482     std::vector<double> fvec;
00483     return fvec;
00484   }
00485 }
00486 
00487 std::vector<std::string> ECalSD::getStringArray(const std::string & str,
00488                                                 const DDsvalues_type & sv) {
00489 
00490 #ifdef DebugLog
00491   LogDebug("EcalSim") << "ECalSD:getStringArray called for " << str;
00492 #endif
00493   DDValue value(str);
00494   if (DDfetch(&sv,value)) {
00495 #ifdef DebugLog
00496     LogDebug("EcalSim") << value;
00497 #endif
00498     const std::vector<std::string> & fvec = value.strings();
00499     return fvec;
00500   } else {
00501     std::vector<std::string> fvec;
00502     return fvec;
00503   }
00504 }