CMS 3D CMS Logo

Public Member Functions | Private Member Functions | Private Attributes

HFShowerPMT Class Reference

#include <HFShowerPMT.h>

List of all members.

Public Member Functions

double getHits (G4Step *aStep)
double getRadius ()
 HFShowerPMT (std::string &name, const DDCompactView &cpv, edm::ParameterSet const &p)
virtual ~HFShowerPMT ()

Private Member Functions

std::vector< double > getDDDArray (const std::string &, const DDsvalues_type &)

Private Attributes

HFCherenkovcherenkov
int indexF
int indexR
double pePerGeV
std::vector< int > pmtFib1
std::vector< int > pmtFib2
std::vector< int > pmtR1
std::vector< int > pmtR2
std::vector< double > rTable

Detailed Description

Definition at line 19 of file HFShowerPMT.h.


Constructor & Destructor Documentation

HFShowerPMT::HFShowerPMT ( std::string &  name,
const DDCompactView cpv,
edm::ParameterSet const &  p 
)

Definition at line 22 of file HFShowerPMT.cc.

References DDFilteredView::addFilter(), cherenkov, DDSpecificsFilter::equals, Exception, DDFilteredView::firstChild(), getDDDArray(), edm::ParameterSet::getParameter(), cuy::ii, getHLTprescales::index, DDFilteredView::mergedSpecifics(), mergeVDriftHistosByStation::name, pePerGeV, pmtFib1, pmtFib2, pmtR1, pmtR2, rTable, DDSpecificsFilter::setCriteria(), and relativeConstraints::value.

                                                    : cherenkov(0) {

  edm::ParameterSet m_HF  = p.getParameter<edm::ParameterSet>("HFShowerPMT");
  pePerGeV                = m_HF.getParameter<double>("PEPerGeVPMT");
  
  G4String attribute = "ReadOutName";
  G4String value     = name;
  DDSpecificsFilter filter0;
  DDValue           ddv0(attribute,value,0);
  filter0.setCriteria(ddv0,DDSpecificsFilter::equals);
  DDFilteredView fv0(cpv);
  fv0.addFilter(filter0);
  if (fv0.firstChild()) {
    DDsvalues_type sv0(fv0.mergedSpecifics());

    //Special Geometry parameters
    rTable   = getDDDArray("rTable",sv0);
    edm::LogInfo("HFShower") << "HFShowerPMT: " << rTable.size() 
                             << " rTable (cm)";
    for (unsigned int ig=0; ig<rTable.size(); ig++)
      edm::LogInfo("HFShower") << "HFShowerPMT: rTable[" << ig << "] = "
                               << rTable[ig]/cm << " cm";
  } else {
    edm::LogError("HFShower") << "HFShowerPMT: cannot get filtered "
                              << " view for " << attribute << " matching "
                              << value;
    throw cms::Exception("Unknown", "HFShowerPMT")
      << "cannot match " << attribute << " to " << name <<"\n";
  }

  attribute = "Volume";
  value     = "HFPMT";
  DDSpecificsFilter filter1;
  DDValue           ddv1(attribute,value,0);
  filter1.setCriteria(ddv1,DDSpecificsFilter::equals);
  DDFilteredView fv1(cpv);
  fv1.addFilter(filter1);
  if (fv1.firstChild()) {
    DDsvalues_type sv1(fv1.mergedSpecifics());
    std::vector<double> neta;
    neta = getDDDArray("indexPMTR",sv1);
    for (unsigned int ii=0; ii<neta.size(); ii++) {
      int index = static_cast<int>(neta[ii]);
      int ir=-1, ifib=-1;
      if (index >= 0) {
        ir   = index/10; ifib = index%10;
      }
      pmtR1.push_back(ir);
      pmtFib1.push_back(ifib);
    }
    neta = getDDDArray("indexPMTL",sv1);
    for (unsigned int ii=0; ii<neta.size(); ii++) {
      int index = static_cast<int>(neta[ii]);
      int ir=-1, ifib=-1;
      if (index >= 0) {
        ir   = index/10; ifib = index%10;
      }
      pmtR2.push_back(ir);
      pmtFib2.push_back(ifib);
    }
    edm::LogInfo("HFShower") << "HFShowerPMT: gets the Index matches for "
                             << neta.size() << " PMTs";
    for (unsigned int ii=0; ii<neta.size(); ii++) 
      edm::LogInfo("HFShower") << "HFShowerPMT: rIndexR[" << ii << "] = "
                               << pmtR1[ii] << " fibreR[" << ii << "] = "
                               << pmtFib1[ii] << " rIndexL[" << ii << "] = "
                               << pmtR2[ii] << " fibreL[" << ii << "] = "
                               << pmtFib2[ii];
  } else {
    edm::LogWarning("HFShower") << "HFShowerPMT: cannot get filtered "
                                << " view for " << attribute << " matching "
                                << value;
  }

  cherenkov = new HFCherenkov(m_HF);
}
HFShowerPMT::~HFShowerPMT ( ) [virtual]

Definition at line 100 of file HFShowerPMT.cc.

References cherenkov.

                          {
  if (cherenkov) delete cherenkov;
}

Member Function Documentation

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

Definition at line 168 of file HFShowerPMT.cc.

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

Referenced by HFShowerPMT().

                                                                        {

#ifdef DebugLog
  LogDebug("HFShower") << "HFShowerPMT:getDDDArray called for " << str;
#endif
  DDValue value(str);
  if (DDfetch(&sv,value)) {
#ifdef DebugLog
    LogDebug("HFShower") << value;
#endif
    const std::vector<double> & fvec = value.doubles();
    int nval = fvec.size();
    if (nval < 2) {
      edm::LogError("HFShower") << "HFShowerPMT: # of " << str 
                                << " bins " << nval << " < 2 ==> illegal";
      throw cms::Exception("Unknown", "HFShowerPMT")
        << "nval < 2 for array " << str << "\n";
    }

    return fvec;
  } else {
    edm::LogError("HFShower") << "HFShowerPMT: cannot get array " << str;
    throw cms::Exception("Unknown", "HFShowerPMT") 
      << "cannot get array " << str << "\n";
  }
}
double HFShowerPMT::getHits ( G4Step *  aStep)

Definition at line 104 of file HFShowerPMT.cc.

References beta, cherenkov, HFCherenkov::computeNPEinPMT(), indexF, indexR, LogDebug, pePerGeV, interactiveExample::photons, pmtFib1, pmtFib2, pmtR1, and pmtR2.

Referenced by HCalSD::getHitPMT().

                                          {

  indexR = indexF = -1;

  G4StepPoint * preStepPoint  = aStep->GetPreStepPoint(); 
  const G4VTouchable* touch   = preStepPoint->GetTouchable();
  int                 boxNo   = touch->GetReplicaNumber(2);
  int                 pmtNo   = touch->GetReplicaNumber(1);
  if (boxNo <= 1) {
    indexR = pmtR1[pmtNo-1];
    indexF = pmtFib1[pmtNo-1];
  } else {
    indexR = pmtR2[pmtNo-1];
    indexF = pmtFib2[pmtNo-1];
  }

#ifdef DebugLog
  double edep = aStep->GetTotalEnergyDeposit();
  LogDebug("HFShower") << "HFShowerPMT: Box " << boxNo << " PMT "
                       << pmtNo << " Mapped Indices " << indexR << ", "
                       << indexF << " Edeposit " << edep/MeV << " MeV; PE "
                       << edep*pePerGeV/GeV;
#endif

  double photons = 0;
  if (indexR >= 0 && indexF > 0) {
    G4Track *aTrack = aStep->GetTrack();
    G4ParticleDefinition *particleDef = aTrack->GetDefinition();
    double stepl = aStep->GetStepLength();
    double beta  = preStepPoint->GetBeta();
    G4ThreeVector pDir = aTrack->GetDynamicParticle()->GetMomentumDirection();
    G4ThreeVector localMom = preStepPoint->GetTouchable()->GetHistory()->
      GetTopTransform().TransformAxis(pDir);
    photons = cherenkov->computeNPEinPMT(particleDef, beta, localMom.x(),
                                         localMom.y(), localMom.z(), stepl);
#ifdef DebugLog
  LogDebug("HFShower") << "HFShowerPMT::getHits: for particle " 
                       << particleDef->GetParticleName() << " Step " << stepl
                       << " Beta " << beta << " Direction " << pDir
                       << " Local " << localMom << " p.e. " << photons;
#endif 

  }
  return photons;
}
double HFShowerPMT::getRadius ( )

Definition at line 150 of file HFShowerPMT.cc.

References indexF, indexR, LogDebug, alignCSCRings::r, and rTable.

Referenced by HCalSD::getHitPMT().

                              {
   
  double r = 0.;
  if (indexR >= 0 && indexR+1 < (int)(rTable.size()))
    r = 0.5*(rTable[indexR]+rTable[indexR+1]);
#ifdef DebugLog
  else
    LogDebug("HFShower") << "HFShowerPMT::getRadius: R " << indexR
                         << " F " << indexF;
#endif
  if (indexF == 2)  r =-r;
#ifdef DebugLog
  LogDebug("HFShower") << "HFShowerPMT: Radius (" << indexR << "/" << indexF 
                       << ") " << r;
#endif
  return r;
}

Member Data Documentation

Definition at line 35 of file HFShowerPMT.h.

Referenced by getHits(), HFShowerPMT(), and ~HFShowerPMT().

int HFShowerPMT::indexF [private]

Definition at line 37 of file HFShowerPMT.h.

Referenced by getHits(), and getRadius().

int HFShowerPMT::indexR [private]

Definition at line 37 of file HFShowerPMT.h.

Referenced by getHits(), and getRadius().

double HFShowerPMT::pePerGeV [private]

Definition at line 36 of file HFShowerPMT.h.

Referenced by getHits(), and HFShowerPMT().

std::vector<int> HFShowerPMT::pmtFib1 [private]

Definition at line 39 of file HFShowerPMT.h.

Referenced by getHits(), and HFShowerPMT().

std::vector<int> HFShowerPMT::pmtFib2 [private]

Definition at line 40 of file HFShowerPMT.h.

Referenced by getHits(), and HFShowerPMT().

std::vector<int> HFShowerPMT::pmtR1 [private]

Definition at line 39 of file HFShowerPMT.h.

Referenced by getHits(), and HFShowerPMT().

std::vector<int> HFShowerPMT::pmtR2 [private]

Definition at line 40 of file HFShowerPMT.h.

Referenced by getHits(), and HFShowerPMT().

std::vector<double> HFShowerPMT::rTable [private]

Definition at line 38 of file HFShowerPMT.h.

Referenced by getRadius(), and HFShowerPMT().