#include <SiTrivialInduceChargeOnStrips.h>
Public Member Functions | |
void | induce (SiChargeCollectionDrifter::collection_type collection_points, const StripGeomDetUnit &det, std::vector< double > &localAmplitudes, size_t &recordMinAffectedStrip, size_t &recordMaxAffectedStrip) |
SiTrivialInduceChargeOnStrips (const edm::ParameterSet &conf, double g) | |
virtual | ~SiTrivialInduceChargeOnStrips () |
Private Member Functions | |
double | chargeDeposited (size_t strip, size_t Nstrips, double amplitude, double chargeSpread, double chargePosition) const |
Static Private Member Functions | |
static unsigned int | indexOf (const std::string &) |
static unsigned int | typeOf (const StripGeomDetUnit &) |
Private Attributes | |
const double | geVperElectron |
const double | Nsigma |
std::vector< std::vector < double > > | signalCoupling |
Static Private Attributes | |
static const int | Ntypes = 14 |
static const std::string | type [] = { "IB1", "IB2","OB1","OB2","W1a","W2a","W3a","W1b","W2b","W3b","W4","W5","W6","W7"} |
Definition at line 7 of file SiTrivialInduceChargeOnStrips.h.
SiTrivialInduceChargeOnStrips::SiTrivialInduceChargeOnStrips | ( | const edm::ParameterSet & | conf, |
double | g | ||
) |
Definition at line 42 of file SiTrivialInduceChargeOnStrips.cc.
References edm::ParameterSet::getParameter(), i, mode, Ntypes, and signalCoupling.
: Nsigma(3.), geVperElectron(g) { std::string mode = conf.getParameter<bool>("APVpeakmode") ? "Peak" : "Dec"; for(int i=0; i<Ntypes; i++) signalCoupling.push_back(conf.getParameter<std::vector<double> >("CouplingConstant"+mode+type[i])); }
virtual SiTrivialInduceChargeOnStrips::~SiTrivialInduceChargeOnStrips | ( | ) | [inline, virtual] |
Definition at line 10 of file SiTrivialInduceChargeOnStrips.h.
{}
double SiTrivialInduceChargeOnStrips::chargeDeposited | ( | size_t | strip, |
size_t | Nstrips, | ||
double | amplitude, | ||
double | chargeSpread, | ||
double | chargePosition | ||
) | const [inline, private] |
Definition at line 89 of file SiTrivialInduceChargeOnStrips.cc.
References geVperElectron.
Referenced by induce().
{ double integralUpToStrip = (strip == 0) ? 0. : ( ROOT::Math::normal_cdf( strip, chargeSpread, chargePosition) ); double integralUpToNext = (strip+1 == Nstrips) ? 1. : ( ROOT::Math::normal_cdf( strip+1, chargeSpread, chargePosition) ); double percentOfSignal = integralUpToNext - integralUpToStrip; return percentOfSignal * amplitude / geVperElectron; }
unsigned int SiTrivialInduceChargeOnStrips::indexOf | ( | const std::string & | t | ) | [inline, static, private] |
void SiTrivialInduceChargeOnStrips::induce | ( | SiChargeCollectionDrifter::collection_type | collection_points, |
const StripGeomDetUnit & | det, | ||
std::vector< double > & | localAmplitudes, | ||
size_t & | recordMinAffectedStrip, | ||
size_t & | recordMaxAffectedStrip | ||
) | [virtual] |
Implements SiInduceChargeOnStrips.
Definition at line 51 of file SiTrivialInduceChargeOnStrips.cc.
References abs, chargeDeposited(), StripTopology::localPitch(), max(), min, Nsigma, StripTopology::nstrips(), signalCoupling, StripGeomDetUnit::specificTopology(), strip(), StripTopology::strip(), and typeOf().
{ std::vector<double>& coupling = signalCoupling.at(typeOf(det)); const StripTopology& topology = dynamic_cast<const StripTopology&>(det.specificTopology()); size_t Nstrips = topology.nstrips(); for (SiChargeCollectionDrifter::collection_type::const_iterator signalpoint = collection_points.begin(); signalpoint != collection_points.end(); signalpoint++ ) { //In strip coordinates: double chargePosition = topology.strip(signalpoint->position()); double chargeSpread = signalpoint->sigma() / topology.localPitch(signalpoint->position()); size_t fromStrip = size_t(std::max( 0, int(std::floor( chargePosition - Nsigma*chargeSpread)))); size_t untilStrip = size_t(std::min( Nstrips, size_t(std::ceil( chargePosition + Nsigma*chargeSpread) ))); for (size_t strip = fromStrip; strip < untilStrip; strip++) { double chargeDepositedOnStrip = chargeDeposited( strip, Nstrips, signalpoint->amplitude(), chargeSpread, chargePosition); size_t affectedFromStrip = size_t(std::max( 0, int(strip - coupling.size() + 1))); size_t affectedUntilStrip = size_t(std::min( Nstrips, strip + coupling.size()) ); for (size_t affectedStrip = affectedFromStrip; affectedStrip < affectedUntilStrip; affectedStrip++) { localAmplitudes.at( affectedStrip ) += chargeDepositedOnStrip * coupling.at(abs( affectedStrip - strip )) ; } if( affectedFromStrip < recordMinAffectedStrip ) recordMinAffectedStrip = affectedFromStrip; if( affectedUntilStrip > recordMaxAffectedStrip ) recordMaxAffectedStrip = affectedUntilStrip; } } return; }
unsigned int SiTrivialInduceChargeOnStrips::typeOf | ( | const StripGeomDetUnit & | det | ) | [inline, static, private] |
Definition at line 30 of file SiTrivialInduceChargeOnStrips.cc.
References Exception, GeomDet::geographicalId(), errorMatrix2Lands_multiChannel::id, indexOf(), TIBDetId::layer(), TIDDetId::ring(), TECDetId::ring(), StripGeomDetUnit::specificType(), GeomDetType::subDetector(), GeomDetEnumerators::TEC, GeomDetEnumerators::TIB, GeomDetEnumerators::TID, and GeomDetEnumerators::TOB.
Referenced by induce().
{ DetId id = det.geographicalId(); switch (det.specificType().subDetector()) { case GeomDetEnumerators::TIB: {return (TIBDetId(id).layer() < 3) ? indexOf("IB1") : indexOf("IB2");} case GeomDetEnumerators::TOB: {return (TOBDetId(id).layer() > 4) ? indexOf("OB1") : indexOf("OB2");} case GeomDetEnumerators::TID: {return indexOf("W1a") -1 + TIDDetId(id).ring();} //fragile: relies on ordering of 'type' case GeomDetEnumerators::TEC: {return indexOf("W1b") -1 + TECDetId(id).ring();} //fragile: relies on ordering of 'type' default: throw cms::Exception("Invalid subdetector") << id(); } }
const double SiTrivialInduceChargeOnStrips::geVperElectron [private] |
Definition at line 30 of file SiTrivialInduceChargeOnStrips.h.
Referenced by chargeDeposited().
const double SiTrivialInduceChargeOnStrips::Nsigma [private] |
Definition at line 29 of file SiTrivialInduceChargeOnStrips.h.
Referenced by induce().
const int SiTrivialInduceChargeOnStrips::Ntypes = 14 [static, private] |
Definition at line 26 of file SiTrivialInduceChargeOnStrips.h.
Referenced by indexOf(), and SiTrivialInduceChargeOnStrips().
std::vector<std::vector<double> > SiTrivialInduceChargeOnStrips::signalCoupling [private] |
Definition at line 27 of file SiTrivialInduceChargeOnStrips.h.
Referenced by induce(), and SiTrivialInduceChargeOnStrips().
const std::string SiTrivialInduceChargeOnStrips::type = { "IB1", "IB2","OB1","OB2","W1a","W2a","W3a","W1b","W2b","W3b","W4","W5","W6","W7"} [static, private] |
Definition at line 25 of file SiTrivialInduceChargeOnStrips.h.
Referenced by indexOf().