CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_7/src/SimTracker/SiStripDigitizer/plugins/SiTrivialInduceChargeOnStrips.cc

Go to the documentation of this file.
00001 #include "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 "DataFormats/Math/interface/approx_erf.h"
00007 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
00008 #include "Geometry/Records/interface/IdealGeometryRecord.h"
00009 #include "FWCore/Utilities/interface/Exception.h"
00010 
00011 #include <algorithm>
00012 #include <iostream>
00013 
00014 namespace {
00015   struct Count {
00016     double ncall=0;
00017     double ndep=0, ndep2=0, maxdep=0;
00018     double nstr=0, nstr2=0;
00019     double ncv=0, nval=0, nval2=0, maxv=0;
00020     double dzero=0;
00021     void dep(double d) { ncall++; ndep+=d; ndep2+=d*d;maxdep=std::max(d,maxdep);}
00022     void str(double d) { nstr+=d; nstr2+=d*d;}
00023     void val(double d) { ncv++; nval+=d; nval2+=d*d; maxv=std::max(d,maxv);}
00024     void zero() { dzero++;}    
00025     ~Count() {
00026 #ifdef SISTRIP_COUNT
00027       std::cout << "deposits " << ncall << " " << maxdep << " " << ndep/ncall << " " << std::sqrt(ndep2*ncall -ndep*ndep)/ncall << std::endl;
00028       std::cout << "zeros " << dzero << std::endl;
00029       std::cout << "strips  " << nstr/ndep << " " << std::sqrt(nstr2*ndep -nstr*nstr)/ndep << std::endl;
00030       std::cout << "vaules  " << ncv << " " << maxv << " " << nval/ncv << " " << std::sqrt(nval2*ncv -nval*nval)/ncv << std::endl;
00031 #endif
00032     }
00033   };
00034   
00035  Count count;
00036 }
00037 
00038 
00039 namespace {
00040   constexpr int Ntypes = 14;
00041   //                                   0     1      2     3     4    5      6     7
00042   const std::string type[Ntypes] = { "IB1", "IB2","OB1","OB2","W1a","W2a","W3a","W1b","W2b","W3b","W4","W5","W6","W7"};
00043   enum { indexOfIB1=0, indexOfIB2=1,  indexOfOB1=2, indexOfOB2=3, indexOfW1a=4, indexOfW1b=7}; 
00044 
00045   inline
00046   std::vector<std::vector<float> >
00047   fillSignalCoupling(const edm::ParameterSet& conf, int nTypes, const std::string* typeArray) {
00048     std::vector<std::vector<float> > signalCoupling;
00049     signalCoupling.reserve(nTypes);
00050     std::string mode = conf.getParameter<bool>("APVpeakmode") ? "Peak" : "Dec";
00051     for(int i=0; i<nTypes; ++i) {
00052       auto dc = conf.getParameter<std::vector<double> >("CouplingConstant"+mode+typeArray[i]);
00053       signalCoupling.emplace_back(dc.begin(),dc.end());
00054     }
00055     return signalCoupling;
00056   }
00057 
00058   inline unsigned int indexOf(const std::string& t) { return std::find( type, type + Ntypes, t) - type;}
00059 
00060 
00061   inline unsigned int typeOf(const StripGeomDetUnit& det, const TrackerTopology *tTopo) {
00062     DetId id = det.geographicalId();
00063     switch (det.specificType().subDetector()) {
00064     case GeomDetEnumerators::TIB: {return (tTopo->tibLayer(id) < 3) ? indexOfIB1 : indexOfIB2;}
00065     case GeomDetEnumerators::TOB: {return (tTopo->tobLayer(id) > 4) ? indexOfOB1 : indexOfOB2;}
00066     case GeomDetEnumerators::TID: {return indexOfW1a -1 + tTopo->tidRing(id);} //fragile: relies on ordering of 'type'
00067     case GeomDetEnumerators::TEC: {return indexOfW1b -1 + tTopo->tecRing(id);} //fragile: relies on ordering of 'type'
00068     default: throw cms::Exception("Invalid subdetector") << id();
00069     }
00070   }
00071 
00072   inline double
00073   chargeDeposited(size_t strip, size_t Nstrips, double amplitude, double chargeSpread, double chargePosition)  {
00074     double integralUpToStrip = (strip == 0)         ? 0. : ( approx_erf(   (strip-chargePosition)/chargeSpread/1.41421356237309515 ) );
00075     double integralUpToNext  = (strip+1 == Nstrips) ? 1. : ( approx_erf(   (strip+1-chargePosition)/chargeSpread/1.41421356237309515 ) );
00076     
00077     return  0.5 * (integralUpToNext - integralUpToStrip) * amplitude;
00078   }
00079   
00080 }
00081 
00082 
00083 SiTrivialInduceChargeOnStrips::
00084 SiTrivialInduceChargeOnStrips(const edm::ParameterSet& conf,double g) 
00085   : signalCoupling(fillSignalCoupling(conf, Ntypes, type)), Nsigma(3.), geVperElectron(g)  {
00086 }
00087 
00088 void 
00089 SiTrivialInduceChargeOnStrips::
00090 induce(const SiChargeCollectionDrifter::collection_type& collection_points, 
00091        const StripGeomDetUnit& det, 
00092        std::vector<float>& localAmplitudes,
00093        size_t& recordMinAffectedStrip, 
00094        size_t& recordMaxAffectedStrip,
00095        const TrackerTopology *tTopo) const {
00096 
00097 
00098    induceVector(collection_points, det, localAmplitudes, recordMinAffectedStrip, recordMaxAffectedStrip, tTopo);
00099 
00100   /*
00101   auto ominA=recordMinAffectedStrip, omaxA=recordMaxAffectedStrip;
00102   std::vector<float> oampl(localAmplitudes);    
00103   induceOriginal(collection_points, det, oampl, ominA, omaxA, tTopo);
00104 
00105   //  std::cout << "orig " << ominA << " " << omaxA << " ";
00106   //for (auto a : oampl) std::cout << a << ",";
00107   //std::cout << std::endl;
00108 
00109   auto minA=recordMinAffectedStrip, maxA=recordMaxAffectedStrip;
00110   std::vector<float> ampl(localAmplitudes);
00111   induceVector(collection_points, det, ampl, minA, maxA, tTopo);
00112 
00113   // std::cout << "vect " << minA << " " << maxA << " ";          
00114   //for (auto a :       ampl) std::cout << a << ",";
00115   //std::cout << std::endl;
00116  
00117   float diff=0;
00118   for (size_t i=0; i!=ampl.size(); ++i) { diff = std::max(diff,ampl[i]>0 ? std::abs(ampl[i]-oampl[i])/ampl[i] : 0);}
00119   if (diff> 1.e-4) {
00120     std::cout << diff << std::endl;
00121     std::cout << "orig " << ominA << " " << omaxA << " ";
00122 //    for (auto a : oampl) std::cout << a << ",";
00123     std::cout << std::endl;
00124     std::cout << "vect " << minA << " " << maxA << " ";
00125 //    for (auto a : ampl) std::cout << a << ",";
00126     std::cout << std::endl;
00127   }
00128 
00129   localAmplitudes.swap(ampl);
00130   recordMinAffectedStrip=minA;
00131   recordMaxAffectedStrip=maxA;
00132   */
00133 
00134 }
00135 
00136 void 
00137 SiTrivialInduceChargeOnStrips::
00138 induceVector(const SiChargeCollectionDrifter::collection_type& collection_points, 
00139                const StripGeomDetUnit& det, 
00140                std::vector<float>& localAmplitudes, 
00141                size_t& recordMinAffectedStrip, 
00142                size_t& recordMaxAffectedStrip,
00143                const TrackerTopology *tTopo) const {
00144 
00145 
00146   auto const & coupling = signalCoupling[typeOf(det,tTopo)];
00147   const StripTopology& topology = dynamic_cast<const StripTopology&>(det.specificTopology());
00148   const int Nstrips =  topology.nstrips();
00149 
00150 
00151   const int NP = collection_points.size();
00152   if(0==NP) return;
00153 
00154   constexpr int MaxN = 512;
00155   // if NP too large split...
00156 
00157   for (int ip=0; ip<NP; ip+=MaxN) {
00158     auto N = std::min(NP-ip,MaxN);
00159 
00160     count.dep(N);
00161     float amplitude[N];
00162     float chargePosition[N];
00163     float chargeSpread[N];
00164     int fromStrip[N];
00165     int nStrip[N];
00166 
00167     // load not vectorize
00168     //In strip coordinates:
00169     for (int i=0; i!=N;++i) {
00170       auto j = ip+i;
00171       if (0==collection_points[j].amplitude()) count.zero();
00172       chargePosition[i]=topology.strip(collection_points[j].position());
00173       chargeSpread[i]= collection_points[j].sigma() / topology.localPitch(collection_points[j].position());
00174       amplitude[i]=0.5f*collection_points[j].amplitude() / geVperElectron;
00175     }
00176     
00177     // this vectorize
00178     for (int i=0; i!=N;++i) {
00179       fromStrip[i]  = std::max( 0,  int(std::floor( chargePosition[i] - Nsigma*chargeSpread[i])) );
00180       nStrip[i] = std::min( Nstrips, int(std::ceil( chargePosition[i] + Nsigma*chargeSpread[i])) ) - fromStrip[i];
00181     }
00182     int tot=0;
00183     for (int i=0; i!=N;++i) tot += nStrip[i];
00184     tot+=N; // add last strip 
00185     count.val(tot);
00186     float value[tot];
00187     
00188     // assign relative position (lower bound of strip) in value;
00189     int kk=0;
00190     for (int i=0; i!=N;++i) {
00191       auto delta = 1.f/(std::sqrt(2.f)*chargeSpread[i]);
00192       auto pos = delta*(float(fromStrip[i])-chargePosition[i]);
00193       for (int j=0;j<=nStrip[i]; ++j)  
00194         value[kk++] = pos+float(j)*delta;  
00195     }
00196     assert(kk==tot);
00197     
00198     // main loop fully vectorized
00199     for (int k=0;k!=tot; ++k)
00200       value[k] = approx_erf(value[k]);
00201     
00202     // saturate 0 & NStrips strip to 0 and 1???
00203     kk=0;
00204     for (int i=0; i!=N;++i) {
00205       if (0 == fromStrip[i])  value[kk]=0;
00206       kk+=nStrip[i];
00207       if (Nstrips == fromStrip[i]+nStrip[i]) value[kk]=1.f;
00208       ++kk;
00209     }
00210     assert(kk==tot);
00211     
00212     // compute integral over strip (lower bound becomes the value)
00213     for (int k=0;k!=tot-1; ++k)
00214       value[k]-=value[k+1];  // this is negative!
00215     
00216     
00217     float charge[Nstrips]; for (int i=0;i!=Nstrips; ++i) charge[i]=0;
00218     kk=0;
00219     for (int i=0; i!=N;++i){ 
00220       for (int j=0;j!=nStrip[i]; ++j)
00221         charge[fromStrip[i]+j]-= amplitude[i]*value[kk++];
00222       ++kk; // skip last "strip"
00223     }
00224     assert(kk==tot);
00225     
00226     
00228     int minA=recordMinAffectedStrip, maxA=recordMaxAffectedStrip;
00229     int sc = coupling.size();
00230     for (int i=0;i!=Nstrips; ++i) {
00231       int strip = i;
00232       if (0==charge[i]) continue;
00233       auto affectedFromStrip  = std::max( 0, strip - sc + 1);
00234       auto affectedUntilStrip = std::min(Nstrips, strip + sc);  
00235       for (auto affectedStrip=affectedFromStrip;  affectedStrip < affectedUntilStrip;  ++affectedStrip)
00236         localAmplitudes[affectedStrip] += charge[i] * coupling[std::abs(affectedStrip - strip)] ;
00237       
00238       if( affectedFromStrip  < minA ) minA = affectedFromStrip;
00239       if( affectedUntilStrip > maxA ) maxA = affectedUntilStrip;
00240     }
00241     recordMinAffectedStrip=minA;
00242     recordMaxAffectedStrip=maxA;
00243   }  // end loop ip
00244 
00245 }
00246 
00247 void 
00248 SiTrivialInduceChargeOnStrips::
00249 induceOriginal(const SiChargeCollectionDrifter::collection_type& collection_points, 
00250                const StripGeomDetUnit& det, 
00251                std::vector<float>& localAmplitudes, 
00252                size_t& recordMinAffectedStrip, 
00253                size_t& recordMaxAffectedStrip,
00254                const TrackerTopology *tTopo) const {
00255 
00256 
00257   auto const & coupling = signalCoupling[typeOf(det,tTopo)];
00258   const StripTopology& topology = dynamic_cast<const StripTopology&>(det.specificTopology());
00259   size_t Nstrips =  topology.nstrips();
00260 
00261   if (!collection_points.empty()) count.dep(collection_points.size());
00262 
00263   for (auto signalpoint = collection_points.begin();  signalpoint != collection_points.end();  signalpoint++ ) {
00264     
00265     //In strip coordinates:
00266     double chargePosition = topology.strip(signalpoint->position());
00267     double chargeSpread = signalpoint->sigma() / topology.localPitch(signalpoint->position());
00268     
00269     size_t fromStrip  = size_t(std::max( 0,          int(std::floor( chargePosition - Nsigma*chargeSpread))));
00270     size_t untilStrip = size_t(std::min( Nstrips, size_t(std::ceil( chargePosition + Nsigma*chargeSpread) )));
00271 
00272     count.str(std::max(0,int(untilStrip)-int(fromStrip)));
00273     for (size_t strip = fromStrip;  strip < untilStrip; strip++) {
00274 
00275       double chargeDepositedOnStrip = chargeDeposited( strip, Nstrips, signalpoint->amplitude() / geVperElectron, chargeSpread, chargePosition);
00276 
00277       size_t affectedFromStrip  = size_t(std::max( 0, int(strip - coupling.size() + 1)));
00278       size_t affectedUntilStrip = size_t(std::min( Nstrips, strip + coupling.size())   );  
00279       for (size_t affectedStrip = affectedFromStrip;  affectedStrip < affectedUntilStrip;  affectedStrip++) {
00280         localAmplitudes.at( affectedStrip ) += chargeDepositedOnStrip * coupling.at(abs( affectedStrip - strip )) ;
00281       }
00282 
00283       if( affectedFromStrip  < recordMinAffectedStrip ) recordMinAffectedStrip = affectedFromStrip;
00284       if( affectedUntilStrip > recordMaxAffectedStrip ) recordMaxAffectedStrip = affectedUntilStrip;
00285     }
00286   }
00287 
00288 }
00289