CMS 3D CMS Logo

Public Member Functions | Private Member Functions | Private Attributes

ECalSD Class Reference

#include <ECalSD.h>

Inheritance diagram for ECalSD:
CaloSD SensitiveCaloDetector Observer< const BeginOfRun * > Observer< const BeginOfEvent * > Observer< const BeginOfTrack * > Observer< const EndOfTrack * > Observer< const EndOfEvent * > SensitiveDetector

List of all members.

Public Member Functions

 ECalSD (G4String, const DDCompactView &, SensitiveDetectorCatalog &, edm::ParameterSet const &p, const SimTrackManager *)
virtual uint16_t getDepth (G4Step *)
virtual double getEnergyDeposit (G4Step *)
virtual uint16_t getRadiationLength (G4Step *)
virtual int getTrackID (G4Track *)
virtual uint32_t setDetUnitId (G4Step *)
void setNumberingScheme (EcalNumberingScheme *)
virtual ~ECalSD ()

Private Member Functions

double crystalLength (G4LogicalVolume *)
double curve_LY (G4Step *)
void getBaseNumber (const G4Step *)
double getBirkL3 (G4Step *)
std::vector< double > getDDDArray (const std::string &, const DDsvalues_type &)
std::vector< std::string > getStringArray (const std::string &, const DDsvalues_type &)
void initMap (G4String, const DDCompactView &)

Private Attributes

EnergyResolutionVsLumi ageing
bool ageingWithSlopeLY
double birk1
double birk2
double birk3
double birkCut
double birkSlope
std::string crystalMat
std::string depth1Name
std::string depth2Name
std::vector< G4LogicalVolume * > noWeight
EcalNumberingSchemenumberingScheme
double slopeLY
bool storeRL
bool storeTrack
EcalBaseNumber theBaseNumber
bool useBirk
bool useBirkL3
std::vector< G4LogicalVolume * > useDepth1
std::vector< G4LogicalVolume * > useDepth2
bool useWeight
std::map< G4LogicalVolume
*, double > 
xtalLMap

Detailed Description

Definition at line 26 of file ECalSD.h.


Constructor & Destructor Documentation

ECalSD::ECalSD ( G4String  name,
const DDCompactView cpv,
SensitiveDetectorCatalog clg,
edm::ParameterSet const &  p,
const SimTrackManager manager 
)

Definition at line 37 of file ECalSD.cc.

References DDFilteredView::addFilter(), ageing, ageingWithSlopeLY, birk1, birk2, birk3, birkCut, birkSlope, crystalMat, depth1Name, depth2Name, DDSpecificsFilter::equals, alcazmumu_cfi::filter, DDFilteredView::firstChild(), g, getDDDArray(), edm::ParameterSet::getParameter(), getStringArray(), edm::ParameterSet::getUntrackedParameter(), initMap(), CaloSD::kmaxIon, CaloSD::kmaxNeutron, CaloSD::kmaxProton, LogDebug, DDFilteredView::mergedSpecifics(), DDSpecificsFilter::setCriteria(), EnergyResolutionVsLumi::setLumies(), setNumberingScheme(), slopeLY, storeRL, storeTrack, AlCaHLTBitMon_QueryRunRegistry::string, CaloSD::suppressHeavy, useBirk, useBirkL3, and useWeight.

                                                                          : 
  CaloSD(name, cpv, clg, p, manager, 
         p.getParameter<edm::ParameterSet>("ECalSD").getParameter<int>("TimeSliceUnit"),
         p.getParameter<edm::ParameterSet>("ECalSD").getParameter<bool>("IgnoreTrackID")), 
  numberingScheme(0){

  //   static SimpleConfigurable<bool>   on1(false,  "ECalSD:UseBirkLaw");
  //   static SimpleConfigurable<double> bk1(0.00463,"ECalSD:BirkC1");
  //   static SimpleConfigurable<double> bk2(-0.03,  "ECalSD:BirkC2");
  //   static SimpleConfigurable<double> bk3(1.0,    "ECalSD:BirkC3");
  // Values from NIM A484 (2002) 239-244: as implemented in Geant3
  //   useBirk          = on1.value();
  //   birk1            = bk1.value()*(g/(MeV*cm2));
  //   birk2            = bk2.value()*(g/(MeV*cm2))*(g/(MeV*cm2));

  edm::ParameterSet m_EC = p.getParameter<edm::ParameterSet>("ECalSD");
  useBirk      = m_EC.getParameter<bool>("UseBirkLaw");
  useBirkL3    = m_EC.getParameter<bool>("BirkL3Parametrization");
  birk1        = m_EC.getParameter<double>("BirkC1")*(g/(MeV*cm2));
  birk2        = m_EC.getParameter<double>("BirkC2");
  birk3        = m_EC.getParameter<double>("BirkC3");
  birkSlope    = m_EC.getParameter<double>("BirkSlope");
  birkCut      = m_EC.getParameter<double>("BirkCut");
  slopeLY      = m_EC.getParameter<double>("SlopeLightYield");
  storeTrack   = m_EC.getParameter<bool>("StoreSecondary");
  crystalMat   = m_EC.getUntrackedParameter<std::string>("XtalMat","E_PbWO4");
  bool isItTB  = m_EC.getUntrackedParameter<bool>("TestBeam", false);
  bool nullNS  = m_EC.getUntrackedParameter<bool>("NullNumbering", false);
  storeRL      = m_EC.getUntrackedParameter<bool>("StoreRadLength", false);
  
  ageingWithSlopeLY   = m_EC.getUntrackedParameter<bool>("AgeingWithSlopeLY", false);
  if(ageingWithSlopeLY) ageing.setLumies(p.getParameter<edm::ParameterSet>("ECalSD").getParameter<double>("DelivLuminosity"),
                                         p.getParameter<edm::ParameterSet>("ECalSD").getParameter<double>("InstLuminosity"));

  //Material list for HB/HE/HO sensitive detectors
  std::string attribute = "ReadOutName";
  DDSpecificsFilter filter;
  DDValue           ddv(attribute,name,0);
  filter.setCriteria(ddv,DDSpecificsFilter::equals);
  DDFilteredView fv(cpv);
  fv.addFilter(filter);
  fv.firstChild();
  DDsvalues_type sv(fv.mergedSpecifics());
  // Use of Weight
  useWeight= true;
  std::vector<double>      tempD = getDDDArray("EnergyWeight",sv);
  if (tempD.size() > 0) { if (tempD[0] < 0.1) useWeight = false; }
  std::vector<std::string> tempS = getStringArray("Depth1Name",sv);
  if (tempS.size() > 0) depth1Name = tempS[0];
  else                  depth1Name = " ";
  tempS = getStringArray("Depth2Name",sv);
  if (tempS.size() > 0) depth2Name = tempS[0];
  else                  depth2Name = " ";

  EcalNumberingScheme* scheme=0;
  if (nullNS)                    scheme = 0;
  else if (name == "EcalHitsEB") scheme = dynamic_cast<EcalNumberingScheme*>(new EcalBarrelNumberingScheme());
  else if (name == "EcalHitsEE") scheme = dynamic_cast<EcalNumberingScheme*>(new EcalEndcapNumberingScheme());
  else if (name == "EcalHitsES") {
    if (isItTB) scheme = dynamic_cast<EcalNumberingScheme*>(new ESTBNumberingScheme());
    else        scheme = dynamic_cast<EcalNumberingScheme*>(new EcalPreshowerNumberingScheme());
    useWeight = false;
  } else {edm::LogWarning("EcalSim") << "ECalSD: ReadoutName not supported\n";}

  if (scheme)  setNumberingScheme(scheme);
#ifdef DebugLog
  LogDebug("EcalSim") << "Constructing a ECalSD  with name " << GetName();
#endif
  if (useWeight) {
    edm::LogInfo("EcalSim")  << "ECalSD:: Use of Birks law is set to      " 
                             << useBirk << "        with three constants kB = "
                             << birk1 << ", C1 = " << birk2 << ", C2 = " 
                             << birk3 <<"\n         Use of L3 parametrization "
                             << useBirkL3 << " with slope " << birkSlope
                             << " and cut off " << birkCut << "\n"
                             << "         Slope for Light yield is set to "
                             << slopeLY;
  } else {
    edm::LogInfo("EcalSim")  << "ECalSD:: energy deposit is not corrected "
                             << " by Birk or light yield curve";
  }

  edm::LogInfo("EcalSim") << "ECalSD:: Suppression Flag " << suppressHeavy
                          << " protons below " << kmaxProton << " MeV,"
                          << " neutrons below " << kmaxNeutron << " MeV and"
                          << " ions below " << kmaxIon << " MeV\n"
                          << "         Depth1 Name = " << depth1Name
                          << " and Depth2 Name = " << depth2Name;

  if (useWeight) initMap(name,cpv);

}
ECalSD::~ECalSD ( ) [virtual]

Definition at line 132 of file ECalSD.cc.

References numberingScheme.

                {
  if (numberingScheme) delete numberingScheme;
}

Member Function Documentation

double ECalSD::crystalLength ( G4LogicalVolume *  lv) [private]

Definition at line 429 of file ECalSD.cc.

References xtalLMap.

Referenced by curve_LY(), and getRadiationLength().

                                                {

  double length= 230.;
  std::map<G4LogicalVolume*,double>::const_iterator ite = xtalLMap.find(lv);
  if (ite != xtalLMap.end()) length = ite->second;
  return length;
}
double ECalSD::curve_LY ( G4Step *  aStep) [private]

Definition at line 388 of file ECalSD.cc.

References ageing, ageingWithSlopeLY, EnergyResolutionVsLumi::calcLightCollectionEfficiencyWeighted(), crystalLength(), CaloSD::currentID, LogDebug, CaloSD::setToLocal(), slopeLY, CaloHitID::unitID(), and histoStyle::weight.

Referenced by getEnergyDeposit().

                                     {

  G4StepPoint*     stepPoint = aStep->GetPreStepPoint();
  G4LogicalVolume* lv        = stepPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();

  double weight = 1.;
  G4ThreeVector  localPoint = setToLocal(stepPoint->GetPosition(),
                                         stepPoint->GetTouchable());

  double crlength = crystalLength(lv);

  if(ageingWithSlopeLY){
    //position along the crystal in mm from 0 to 230 (in EB)
    double depth = 0.5 * crlength + localPoint.z();

    if (depth >= -0.1 || depth <= crlength+0.1)
      weight = ageing.calcLightCollectionEfficiencyWeighted(currentID.unitID(), depth/crlength);
  }
  else{
    double dapd = 0.5 * crlength - localPoint.z();
    if (dapd >= -0.1 || dapd <= crlength+0.1) {
      if (dapd <= 100.)
        weight = 1.0 + slopeLY - dapd * 0.01 * slopeLY;
    } else {
      edm::LogWarning("EcalSim") << "ECalSD: light coll curve : wrong distance "
                                 << "to APD " << dapd << " crlength = " 
                                 << crlength <<" crystal name = " <<lv->GetName()
                                 << " z of localPoint = " << localPoint.z() 
                                 << " take weight = " << weight;
    }
  }
#ifdef DebugLog
  LogDebug("EcalSim") << "ECalSD, light coll curve : " << dapd 
                      << " crlength = " << crlength
                      << " crystal name = " << lv->GetName()
                      << " z of localPoint = " << localPoint.z() 
                      << " take weight = " << weight;
#endif
  return weight;
}
void ECalSD::getBaseNumber ( const G4Step *  aStep) [private]

Definition at line 437 of file ECalSD.cc.

References EcalBaseNumber::addLevel(), EcalBaseNumber::getCapacity(), cuy::ii, LogDebug, EcalBaseNumber::reset(), EcalBaseNumber::setSize(), and theBaseNumber.

Referenced by setDetUnitId().

                                              {

  theBaseNumber.reset();
  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
  int theSize = touch->GetHistoryDepth()+1;
  if ( theBaseNumber.getCapacity() < theSize ) theBaseNumber.setSize(theSize);
  //Get name and copy numbers
  if ( theSize > 1 ) {
    for (int ii = 0; ii < theSize ; ii++) {
      theBaseNumber.addLevel(touch->GetVolume(ii)->GetName(),touch->GetReplicaNumber(ii));
#ifdef DebugLog
      LogDebug("EcalSim") << "ECalSD::getBaseNumber(): Adding level " << ii
                          << ": " << touch->GetVolume(ii)->GetName() << "["
                          << touch->GetReplicaNumber(ii) << "]";
#endif
    }
  }

}
double ECalSD::getBirkL3 ( G4Step *  aStep) [private]

Definition at line 457 of file ECalSD.cc.

References birk1, birkCut, birkSlope, DeDxDiscriminatorTools::charge(), create_public_lumi_plots::log, LogDebug, and histoStyle::weight.

Referenced by getEnergyDeposit().

                                      {

  double weight = 1.;
  double charge = aStep->GetPreStepPoint()->GetCharge();

  if (charge != 0. && aStep->GetStepLength() > 0) {
    G4Material* mat = aStep->GetPreStepPoint()->GetMaterial();
    double density = mat->GetDensity();
    double dedx    = aStep->GetTotalEnergyDeposit()/aStep->GetStepLength();
    double rkb     = birk1/density;
    if (dedx > 0) {
      weight         = 1. - birkSlope*log(rkb*dedx);
      if (weight < birkCut) weight = birkCut;
      else if (weight > 1.) weight = 1.;
    }
#ifdef DebugLog
    LogDebug("EcalSim") << "ECalSD::getBirkL3 in " << mat->GetName()
                        << " Charge " << charge << " dE/dx " << dedx
                        << " Birk Const " << rkb << " Weight = " << weight 
                        << " dE " << aStep->GetTotalEnergyDeposit();
#endif
  }
  return weight;

}
std::vector< double > ECalSD::getDDDArray ( const std::string &  str,
const DDsvalues_type sv 
) [private]

Definition at line 483 of file ECalSD.cc.

References DDfetch(), DDValue::doubles(), LogDebug, and relativeConstraints::value.

Referenced by ECalSD().

                                                                   {

#ifdef DebugLog
  LogDebug("EcalSim") << "ECalSD:getDDDArray called for " << str;
#endif
  DDValue value(str);
  if (DDfetch(&sv,value)) {
#ifdef DebugLog
    LogDebug("EcalSim") << value;
#endif
    const std::vector<double> & fvec = value.doubles();
    return fvec;
  } else {
    std::vector<double> fvec;
    return fvec;
  }
}
uint16_t ECalSD::getDepth ( G4Step *  aStep) [virtual]

Reimplemented from CaloSD.

Definition at line 224 of file ECalSD.cc.

References any, getRadiationLength(), LogDebug, run_regression::ret, storeRL, useDepth1, and useDepth2.

                                        {
  G4LogicalVolume* lv   = aStep->GetPreStepPoint()->GetTouchable()->GetVolume(0)->GetLogicalVolume();
  uint16_t ret = 0;
  if (any(useDepth1,lv))      ret = 1;
  else if (any(useDepth2,lv)) ret = 2;
  else if (storeRL) ret = getRadiationLength(aStep);
#ifdef DebugLog
  LogDebug("EcalSim") << "Volume " << lv->GetName() << " Depth " << ret;
#endif
  return ret;
}
double ECalSD::getEnergyDeposit ( G4Step *  aStep) [virtual]

Reimplemented from CaloSD.

Definition at line 136 of file ECalSD.cc.

References any, birk1, birk2, birk3, curve_LY(), CaloSD::getAttenuation(), getBirkL3(), CaloSD::getResponseWt(), TrackInformation::isPrimary(), ke, CaloSD::kmaxIon, CaloSD::kmaxNeutron, CaloSD::kmaxProton, LogDebug, noWeight, NULL, CaloSD::preStepPoint, CaloSD::suppressHeavy, CaloSD::theTrack, useBirk, useBirkL3, useWeight, and histoStyle::weight.

                                              {
  
  if (aStep == NULL) {
    return 0;
  } else {
    preStepPoint        = aStep->GetPreStepPoint();
    G4Track* theTrack   = aStep->GetTrack();
    double wt2  = theTrack->GetWeight();
    G4String nameVolume = preStepPoint->GetPhysicalVolume()->GetName();

    // take into account light collection curve for crystals
    double weight = 1.;
    if (suppressHeavy) {
      TrackInformation * trkInfo = (TrackInformation *)(theTrack->GetUserInformation());
      if (trkInfo) {
        int pdg = theTrack->GetDefinition()->GetPDGEncoding();
        if (!(trkInfo->isPrimary())) { // Only secondary particles
          double ke = theTrack->GetKineticEnergy()/MeV;
          if (((pdg/1000000000 == 1 && ((pdg/10000)%100) > 0 &&
                ((pdg/10)%100) > 0)) && (ke<kmaxIon)) weight = 0;
          if ((pdg == 2212) && (ke < kmaxProton))     weight = 0;
          if ((pdg == 2112) && (ke < kmaxNeutron))    weight = 0;
#ifdef DebugLog
          if (weight == 0) 
            LogDebug("EcalSim") << "Ignore Track " << theTrack->GetTrackID()
                                << " Type " << theTrack->GetDefinition()->GetParticleName()
                                << " Kinetic Energy " << ke << " MeV";
#endif
        }
      }
    }
    G4LogicalVolume* lv   = aStep->GetPreStepPoint()->GetTouchable()->GetVolume(0)->GetLogicalVolume();
    if (useWeight && !any(noWeight,lv)) {
      weight *= curve_LY(aStep);
      if (useBirk) {
        if (useBirkL3) weight *= getBirkL3(aStep);
        else           weight *= getAttenuation(aStep, birk1, birk2, birk3);
      }
    }
    double wt1  = getResponseWt(theTrack);
    double edep = aStep->GetTotalEnergyDeposit()*weight*wt1;
    /*
    if(wt2 != 1.0) { 
      std::cout << "ECalSD:: " << nameVolume
                <<" LightWeight= " <<weight << " wt1= " <<wt1
                << "  wt2= " << wt2 << "  "
                << " Weighted Energy Deposit " << edep/MeV << " MeV" 
                << std::endl;
      std::cout << theTrack->GetDefinition()->GetParticleName()
                << " " << theTrack->GetKineticEnergy()
                << " Id=" << theTrack->GetTrackID()
                << " IdP=" << theTrack->GetParentID();
      const G4VProcess* pr = theTrack->GetCreatorProcess();
      if(pr) std::cout << " from  " << pr->GetProcessName();
      std::cout << std::endl;
    }
    */
    if(wt2 > 0.0) { edep *= wt2; }
#ifdef DebugLog
    LogDebug("EcalSim") << "ECalSD:: " << nameVolume
                        <<" Light Collection Efficiency " <<weight << ":" <<wt1
                        << " Weighted Energy Deposit " << edep/MeV << " MeV";
#endif
    return edep;
  } 
}
uint16_t ECalSD::getRadiationLength ( G4Step *  aStep) [virtual]

Definition at line 236 of file ECalSD.cc.

References crystalLength(), NULL, CaloSD::setToLocal(), and useWeight.

Referenced by getDepth().

                                                  {
  
  uint16_t thisX0 = 0;
  if (aStep != NULL) {
    G4StepPoint* hitPoint = aStep->GetPreStepPoint();
    G4LogicalVolume* lv   = hitPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
    
    if (useWeight) {
      G4ThreeVector  localPoint = setToLocal(hitPoint->GetPosition(),
                                             hitPoint->GetTouchable());
      double crlength = crystalLength(lv);
      double radl     = hitPoint->GetMaterial()->GetRadlen();
      double detz     = (float)(0.5*crlength + localPoint.z());
      thisX0 = (uint16_t)floor(detz/radl);   
    } 
  }
  return thisX0;
}
std::vector< std::string > ECalSD::getStringArray ( const std::string &  str,
const DDsvalues_type sv 
) [private]

Definition at line 502 of file ECalSD.cc.

References DDfetch(), LogDebug, DDValue::strings(), and relativeConstraints::value.

Referenced by ECalSD().

                                                                           {

#ifdef DebugLog
  LogDebug("EcalSim") << "ECalSD:getStringArray called for " << str;
#endif
  DDValue value(str);
  if (DDfetch(&sv,value)) {
#ifdef DebugLog
    LogDebug("EcalSim") << value;
#endif
    const std::vector<std::string> & fvec = value.strings();
    return fvec;
  } else {
    std::vector<std::string> fvec;
    return fvec;
  }
}
int ECalSD::getTrackID ( G4Track *  aTrack) [virtual]

Reimplemented from CaloSD.

Definition at line 203 of file ECalSD.cc.

References any, CaloSD::forceSave, CaloSD::preStepPoint, storeTrack, useDepth1, and useDepth2.

                                      {

  int  primaryID(0);
  bool flag(false);
  if (storeTrack) {
    G4LogicalVolume* lv  = preStepPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
    if (any(useDepth1,lv)) {
      flag = true;
    } else if (any(useDepth2,lv)) {
      flag = true;
    }
  }
  if (flag) {
    forceSave = true;
    primaryID = aTrack->GetTrackID();
  } else {
    primaryID = CaloSD::getTrackID(aTrack);
  }
  return primaryID;
}
void ECalSD::initMap ( G4String  sd,
const DDCompactView cpv 
) [private]

Definition at line 274 of file ECalSD.cc.

References DDFilteredView::addFilter(), any, crystalMat, ddtrap, depth1Name, depth2Name, DDSpecificsFilter::equals, alcazmumu_cfi::filter, DDFilteredView::firstChild(), i, LogDebug, DDFilteredView::logicalPart(), DDLogicalPart::material(), SensitiveDetector::name, DDName::name(), DDBase< N, C >::name(), DDFilteredView::next(), noWeight, DDSolid::parameters(), DDSpecificsFilter::setCriteria(), DDSolid::shape(), DDLogicalPart::solid(), AlCaHLTBitMon_QueryRunRegistry::string, useDepth1, useDepth2, and xtalLMap.

Referenced by ECalSD().

                                                           {

  G4String attribute = "ReadOutName";
  DDSpecificsFilter filter;
  DDValue           ddv(attribute,sd,0);
  filter.setCriteria(ddv,DDSpecificsFilter::equals);
  DDFilteredView fv(cpv);
  fv.addFilter(filter);
  fv.firstChild();

  std::vector<G4LogicalVolume*> lvused;
  const G4LogicalVolumeStore *  lvs = G4LogicalVolumeStore::GetInstance();
  std::vector<G4LogicalVolume *>::const_iterator lvcite;
  std::map<std::string, G4LogicalVolume *> nameMap;
  for (auto lvi = lvs->begin(), lve = lvs->end(); lvi != lve; ++lvi)
    nameMap.insert(std::make_pair((*lvi)->GetName(), *lvi));

  bool dodet=true;
  while (dodet) {
    const std::string &matname = fv.logicalPart().material().name().name();
    const std::string &lvname = fv.logicalPart().name().name();
    G4LogicalVolume* lv = nameMap[lvname];
    if (depth1Name != " ") {
      if (strncmp(lvname.c_str(), depth1Name.c_str(), 4) == 0) {
        if (!any(useDepth1, lv)) {
          useDepth1.push_back(lv);
#ifdef DebugLog
          LogDebug("EcalSim") << "ECalSD::initMap Logical Volume " << lvname
                              <<" in Depth 1 volume list";
#endif
        }
        G4LogicalVolume* lvr = nameMap[lvname + "_refl"];
        if (lvr != 0 && !any(useDepth1, lvr)) {
          useDepth1.push_back(lvr);
#ifdef DebugLog
          LogDebug("EcalSim") << "ECalSD::initMap Logical Volume " << lvname << "_refl"
                              <<" in Depth 1 volume list";
#endif
        }
      }
    }
    if (depth2Name != " ") {
      if (strncmp(lvname.c_str(), depth2Name.c_str(), 4) == 0) {
        if (!any(useDepth2, lv)) {
          useDepth2.push_back(lv);
#ifdef DebugLog
          LogDebug("EcalSim") << "ECalSD::initMap Logical Volume " << lvname
                              <<" in Depth 2 volume list";
#endif
        }
        G4LogicalVolume* lvr = nameMap[lvname + "_refl"];
        if (lvr != 0 && !any(useDepth2,lvr)) {
          useDepth2.push_back(lvr);
#ifdef DebugLog
          LogDebug("EcalSim") << "ECalSD::initMap Logical Volume " << lvname << "_refl"
                              <<" in Depth 2 volume list";
#endif
        }
      }
    }
    if (lv != 0) {
      if (crystalMat.size() == matname.size() && !strcmp(crystalMat.c_str(), matname.c_str())) {
        if (!any(lvused,lv)) {
          lvused.push_back(lv);
          const DDSolid & sol  = fv.logicalPart().solid();
          const std::vector<double> & paras = sol.parameters();
#ifdef DebugLog
          LogDebug("EcalSim") << "ECalSD::initMap (for " << sd << "): Solid " 
                              << lvname << " Shape " << sol.shape() 
                              << " Parameter 0 = "<< paras[0] 
                              << " Logical Volume " << lv;
#endif
          if (sol.shape() == ddtrap) {
            double dz = 2*paras[0];
            xtalLMap.insert(std::pair<G4LogicalVolume*,double>(lv,dz));
            lv = nameMap[lvname + "_refl"];
            if (lv != 0)
              xtalLMap.insert(std::pair<G4LogicalVolume*,double>(lv,dz));
          }
        }
      } else {
        if (!any(noWeight,lv)) {
          noWeight.push_back(lv);
#ifdef DebugLog
          LogDebug("EcalSim") << "ECalSD::initMap Logical Volume " << lvname
                              << " Material " << matname <<" in noWeight list";
#endif
        }
        lv = nameMap[lvname];
        if (lv != 0 && !any(noWeight,lv)) {
          noWeight.push_back(lv);
#ifdef DebugLog
          LogDebug("EcalSim") << "ECalSD::initMap Logical Volume " << lvname
                              << " Material " << matname <<" in noWeight list";
#endif
        }
      }
    }
    dodet = fv.next();
  }
#ifdef DebugLog
  LogDebug("EcalSim") << "ECalSD: Length Table for " << attribute << " = " 
                      << sd << ":";   
  std::map<G4LogicalVolume*,double>::const_iterator ite = xtalLMap.begin();
  int i=0;
  for (; ite != xtalLMap.end(); ite++, i++) {
    G4String name = "Unknown";
    if (ite->first != 0) name = (ite->first)->GetName();
    LogDebug("EcalSim") << " " << i << " " << ite->first << " " << name 
                        << " L = " << ite->second;
  }
#endif
}
uint32_t ECalSD::setDetUnitId ( G4Step *  aStep) [virtual]

Implements CaloSD.

Definition at line 255 of file ECalSD.cc.

References getBaseNumber(), EcalNumberingScheme::getUnitID(), numberingScheme, and theBaseNumber.

                                            { 
  if (numberingScheme == 0) {
    return EBDetId(1,1)();
  } else {
    getBaseNumber(aStep);
    return numberingScheme->getUnitID(theBaseNumber);
  }
}
void ECalSD::setNumberingScheme ( EcalNumberingScheme scheme)

Definition at line 264 of file ECalSD.cc.

References numberingScheme.

Referenced by ECalSD(), and HcalTB04Analysis::update().

                                                           {
  if (scheme != 0) {
    edm::LogInfo("EcalSim") << "EcalSD: updates numbering scheme for " 
                            << GetName() << "\n";
    if (numberingScheme) delete numberingScheme;
    numberingScheme = scheme;
  }
}

Member Data Documentation

Definition at line 61 of file ECalSD.h.

Referenced by curve_LY(), and ECalSD().

bool ECalSD::ageingWithSlopeLY [private]

Definition at line 62 of file ECalSD.h.

Referenced by curve_LY(), and ECalSD().

double ECalSD::birk1 [private]

Definition at line 55 of file ECalSD.h.

Referenced by ECalSD(), getBirkL3(), and getEnergyDeposit().

double ECalSD::birk2 [private]

Definition at line 55 of file ECalSD.h.

Referenced by ECalSD(), and getEnergyDeposit().

double ECalSD::birk3 [private]

Definition at line 55 of file ECalSD.h.

Referenced by ECalSD(), and getEnergyDeposit().

double ECalSD::birkCut [private]

Definition at line 55 of file ECalSD.h.

Referenced by ECalSD(), and getBirkL3().

double ECalSD::birkSlope [private]

Definition at line 55 of file ECalSD.h.

Referenced by ECalSD(), and getBirkL3().

std::string ECalSD::crystalMat [private]

Definition at line 57 of file ECalSD.h.

Referenced by ECalSD(), and initMap().

std::string ECalSD::depth1Name [private]

Definition at line 57 of file ECalSD.h.

Referenced by ECalSD(), and initMap().

std::string ECalSD::depth2Name [private]

Definition at line 57 of file ECalSD.h.

Referenced by ECalSD(), and initMap().

std::vector<G4LogicalVolume*> ECalSD::noWeight [private]

Definition at line 59 of file ECalSD.h.

Referenced by getEnergyDeposit(), and initMap().

Definition at line 52 of file ECalSD.h.

Referenced by setDetUnitId(), setNumberingScheme(), and ~ECalSD().

double ECalSD::slopeLY [private]

Definition at line 56 of file ECalSD.h.

Referenced by curve_LY(), and ECalSD().

bool ECalSD::storeRL [private]

Definition at line 53 of file ECalSD.h.

Referenced by ECalSD(), and getDepth().

bool ECalSD::storeTrack [private]

Definition at line 53 of file ECalSD.h.

Referenced by ECalSD(), and getTrackID().

Definition at line 60 of file ECalSD.h.

Referenced by getBaseNumber(), and setDetUnitId().

bool ECalSD::useBirk [private]

Definition at line 54 of file ECalSD.h.

Referenced by ECalSD(), and getEnergyDeposit().

bool ECalSD::useBirkL3 [private]

Definition at line 54 of file ECalSD.h.

Referenced by ECalSD(), and getEnergyDeposit().

std::vector<G4LogicalVolume*> ECalSD::useDepth1 [private]

Definition at line 59 of file ECalSD.h.

Referenced by getDepth(), getTrackID(), and initMap().

std::vector<G4LogicalVolume*> ECalSD::useDepth2 [private]

Definition at line 59 of file ECalSD.h.

Referenced by getDepth(), getTrackID(), and initMap().

bool ECalSD::useWeight [private]

Definition at line 53 of file ECalSD.h.

Referenced by ECalSD(), getEnergyDeposit(), and getRadiationLength().

std::map<G4LogicalVolume*,double> ECalSD::xtalLMap [private]

Definition at line 58 of file ECalSD.h.

Referenced by crystalLength(), and initMap().