CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_2_9_HLT1_bphpatch4/src/SimTracker/SiStripDigitizer/src/SiTrivialInduceChargeOnStrips.cc

Go to the documentation of this file.
00001 #include "SimTracker/SiStripDigitizer/interface/SiTrivialInduceChargeOnStrips.h"
00002 
00003 #include "Geometry/TrackerGeometryBuilder/interface/StripGeomDetUnit.h"
00004 #include "Geometry/CommonTopologies/interface/StripTopology.h"
00005 #include "Geometry/TrackerGeometryBuilder/interface/StripGeomDetType.h"
00006 #include <Math/ProbFuncMathCore.h>
00007 #include "DataFormats/SiStripDetId/interface/TECDetId.h"
00008 #include "DataFormats/SiStripDetId/interface/TIBDetId.h"
00009 #include "DataFormats/SiStripDetId/interface/TIDDetId.h"
00010 #include "DataFormats/SiStripDetId/interface/TOBDetId.h"
00011 #include "FWCore/Utilities/interface/Exception.h"
00012 
00013 #include <algorithm>
00014 #include <iostream>
00015 
00016 const int
00017 SiTrivialInduceChargeOnStrips::
00018 Ntypes = 14;
00019 
00020 const  std::string 
00021 SiTrivialInduceChargeOnStrips::
00022 type[Ntypes] = { "IB1", "IB2","OB1","OB2","W1a","W2a","W3a","W1b","W2b","W3b","W4","W5","W6","W7"};
00023 
00024 inline unsigned int 
00025 SiTrivialInduceChargeOnStrips::
00026 indexOf(const std::string& t) { return std::find( type, type + Ntypes, t) - type;}
00027 
00028 inline unsigned int
00029 SiTrivialInduceChargeOnStrips::
00030 typeOf(const StripGeomDetUnit& det) {
00031   DetId id = det.geographicalId();
00032   switch (det.specificType().subDetector()) {
00033   case GeomDetEnumerators::TIB: {return (TIBDetId(id).layer() < 3) ? indexOf("IB1") : indexOf("IB2");}
00034   case GeomDetEnumerators::TOB: {return (TOBDetId(id).layer() > 4) ? indexOf("OB1") : indexOf("OB2");}
00035   case GeomDetEnumerators::TID: {return indexOf("W1a") -1 + TIDDetId(id).ring();} //fragile: relies on ordering of 'type'
00036   case GeomDetEnumerators::TEC: {return indexOf("W1b") -1 + TECDetId(id).ring();} //fragile: relies on ordering of 'type'
00037   default: throw cms::Exception("Invalid subdetector") << id();
00038   }
00039 }
00040 
00041 SiTrivialInduceChargeOnStrips::
00042 SiTrivialInduceChargeOnStrips(const edm::ParameterSet& conf,double g) 
00043   : Nsigma(3.), geVperElectron(g)  {
00044   std::string mode = conf.getParameter<bool>("APVpeakmode") ? "Peak" : "Dec";
00045   for(int i=0; i<Ntypes; i++)
00046     signalCoupling.push_back(conf.getParameter<std::vector<double> >("CouplingConstant"+mode+type[i]));
00047 }
00048 
00049 void 
00050 SiTrivialInduceChargeOnStrips::
00051 induce(SiChargeCollectionDrifter::collection_type collection_points, 
00052        const StripGeomDetUnit& det, 
00053        std::vector<double>& localAmplitudes, 
00054        size_t& recordMinAffectedStrip, 
00055        size_t& recordMaxAffectedStrip) {
00056 
00057   std::vector<double>& coupling = signalCoupling.at(typeOf(det));
00058   const StripTopology& topology = dynamic_cast<const StripTopology&>(det.specificTopology());
00059   size_t Nstrips =  topology.nstrips();
00060 
00061   for (SiChargeCollectionDrifter::collection_type::const_iterator 
00062          signalpoint = collection_points.begin();  signalpoint != collection_points.end();  signalpoint++ ) {
00063     
00064     //In strip coordinates:
00065     double chargePosition = topology.strip(signalpoint->position());
00066     double chargeSpread = signalpoint->sigma() / topology.localPitch(signalpoint->position());
00067     
00068     size_t fromStrip  = size_t(std::max( 0,          int(std::floor( chargePosition - Nsigma*chargeSpread))));
00069     size_t untilStrip = size_t(std::min( Nstrips, size_t(std::ceil( chargePosition + Nsigma*chargeSpread) )));
00070     for (size_t strip = fromStrip;  strip < untilStrip; strip++) {
00071 
00072       double chargeDepositedOnStrip = chargeDeposited( strip, Nstrips, signalpoint->amplitude(), chargeSpread, chargePosition);
00073 
00074       size_t affectedFromStrip  = size_t(std::max( 0, int(strip - coupling.size() + 1)));
00075       size_t affectedUntilStrip = size_t(std::min( Nstrips, strip + coupling.size())   );  
00076       for (size_t affectedStrip = affectedFromStrip;  affectedStrip < affectedUntilStrip;  affectedStrip++) {
00077         localAmplitudes.at( affectedStrip ) += chargeDepositedOnStrip * coupling.at(abs( affectedStrip - strip )) ;
00078       }
00079 
00080       if( affectedFromStrip  < recordMinAffectedStrip ) recordMinAffectedStrip = affectedFromStrip;
00081       if( affectedUntilStrip > recordMaxAffectedStrip ) recordMaxAffectedStrip = affectedUntilStrip;
00082     }
00083   }
00084   return;
00085 }
00086 
00087 inline double
00088 SiTrivialInduceChargeOnStrips::
00089 chargeDeposited(size_t strip, size_t Nstrips, double amplitude, double chargeSpread, double chargePosition) const {
00090   double integralUpToStrip = (strip == 0)         ? 0. : ( ROOT::Math::normal_cdf(   strip, chargeSpread, chargePosition) );
00091   double integralUpToNext  = (strip+1 == Nstrips) ? 1. : ( ROOT::Math::normal_cdf( strip+1, chargeSpread, chargePosition) );
00092   double percentOfSignal = integralUpToNext - integralUpToStrip;
00093   
00094   return percentOfSignal * amplitude / geVperElectron;
00095 }