00001
00002
00003
00004 #include <vector>
00005 #include <algorithm>
00006 #include <iostream>
00007 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00008 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00009 #include "SimTracker/SiStripDigitizer/interface/SiStripDigitizerAlgorithm.h"
00010 #include "SimDataFormats/TrackingHit/interface/PSimHitContainer.h"
00011 #include "DataFormats/Common/interface/DetSetVector.h"
00012 #include "DataFormats/GeometrySurface/interface/BoundSurface.h"
00013 #include "CLHEP/Random/RandFlat.h"
00014
00015 #define CBOLTZ (1.38E-23)
00016 #define e_SI (1.6E-19)
00017
00018 SiStripDigitizerAlgorithm::SiStripDigitizerAlgorithm(const edm::ParameterSet& conf, CLHEP::HepRandomEngine& eng):
00019 conf_(conf),rndEngine(eng){
00020 theThreshold = conf_.getParameter<double>("NoiseSigmaThreshold");
00021 theFedAlgo = conf_.getParameter<int>("FedAlgorithm");
00022 peakMode = conf_.getParameter<bool>("APVpeakmode");
00023 theElectronPerADC = conf_.getParameter<double>( peakMode ? "electronPerAdcPeak" : "electronPerAdcDec" );
00024 noise = conf_.getParameter<bool>("Noise");
00025 zeroSuppression = conf_.getParameter<bool>("ZeroSuppression");
00026 theTOFCutForPeak = conf_.getParameter<double>("TOFCutForPeak");
00027 theTOFCutForDeconvolution = conf_.getParameter<double>("TOFCutForDeconvolution");
00028 cosmicShift = conf_.getUntrackedParameter<double>("CosmicDelayShift");
00029 inefficiency = conf_.getParameter<double>("Inefficiency");
00030 RealPedestals = conf_.getParameter<bool>("RealPedestals");
00031 SingleStripNoise = conf_.getParameter<bool>("SingleStripNoise");
00032 CommonModeNoise = conf_.getParameter<bool>("CommonModeNoise");
00033 BaselineShift = conf_.getParameter<bool>("BaselineShift");
00034 APVSaturationFromHIP = conf_.getParameter<bool>("APVSaturationFromHIP");
00035 APVSaturationProb = conf_.getParameter<double>("APVSaturationProb");
00036 cmnRMStib = conf_.getParameter<double>("cmnRMStib");
00037 cmnRMStob = conf_.getParameter<double>("cmnRMStob");
00038 cmnRMStid = conf_.getParameter<double>("cmnRMStid");
00039 cmnRMStec = conf_.getParameter<double>("cmnRMStec");
00040 pedOffset = (unsigned int)conf_.getParameter<double>("PedestalsOffset");
00041 if (peakMode) {
00042 tofCut=theTOFCutForPeak;
00043 LogDebug("StripDigiInfo")<<"APVs running in peak mode (poor time resolution)";
00044 } else {
00045 tofCut=theTOFCutForDeconvolution;
00046 LogDebug("StripDigiInfo")<<"APVs running in deconvolution mode (good time resolution)";
00047 };
00048
00049 theSiHitDigitizer = new SiHitDigitizer(conf_,rndEngine);
00050 theSiPileUpSignals = new SiPileUpSignals();
00051 theSiNoiseAdder = new SiGaussianTailNoiseAdder(theThreshold,rndEngine);
00052 theSiDigitalConverter = new SiTrivialDigitalConverter(theElectronPerADC);
00053 theSiZeroSuppress = new SiStripFedZeroSuppression(theFedAlgo);
00054 theFlatDistribution = new CLHEP::RandFlat(rndEngine, 0., 1.);
00055
00056
00057 }
00058
00059 SiStripDigitizerAlgorithm::~SiStripDigitizerAlgorithm(){
00060 delete theSiHitDigitizer;
00061 delete theSiPileUpSignals;
00062 delete theSiNoiseAdder;
00063 delete theSiDigitalConverter;
00064 delete theSiZeroSuppress;
00065 delete theFlatDistribution;
00066
00067
00068 }
00069
00070
00071
00072
00073 void SiStripDigitizerAlgorithm::run(edm::DetSet<SiStripDigi>& outdigi,
00074 edm::DetSet<SiStripRawDigi>& outrawdigi,
00075 const std::vector<std::pair<const PSimHit*, int > > &input,
00076 StripGeomDetUnit *det,
00077 GlobalVector bfield,float langle,
00078 edm::ESHandle<SiStripGain> & gainHandle,
00079 edm::ESHandle<SiStripThreshold> & thresholdHandle,
00080 edm::ESHandle<SiStripNoises> & noiseHandle,
00081 edm::ESHandle<SiStripPedestals> & pedestalHandle,
00082 edm::ESHandle<SiStripBadStrip> & deadChannelHandle
00083 ) {
00084 theSiPileUpSignals->reset();
00085 unsigned int detID = det->geographicalId().rawId();
00086 SiStripNoises::Range detNoiseRange = noiseHandle->getRange(detID);
00087 SiStripApvGain::Range detGainRange = gainHandle->getRange(detID);
00088 SiStripPedestals::Range detPedestalRange = pedestalHandle->getRange(detID);
00089 SiStripBadStrip::Range detBadStripRange = deadChannelHandle->getRange(detID);
00090 numStrips = (det->specificTopology()).nstrips();
00091
00092
00093 std::vector<bool> badChannels;
00094 badChannels.clear();
00095 badChannels.insert(badChannels.begin(),numStrips,false);
00096 SiStripBadStrip::data fs;
00097 for(SiStripBadStrip::ContainerIterator it=detBadStripRange.first;it!=detBadStripRange.second;++it){
00098 fs=deadChannelHandle->decode(*it);
00099 for(int strip = fs.firstStrip; strip <fs.firstStrip+fs.range; ++strip )badChannels[strip] = true;
00100 }
00101
00102
00103
00104 locAmpl.clear();
00105 detAmpl.clear();
00106 locAmpl.insert(locAmpl.begin(),numStrips,0.);
00107
00108 detAmpl.insert(detAmpl.begin(),numStrips,0.);
00109
00110 firstChannelWithSignal = numStrips;
00111 lastChannelWithSignal = 0;
00112
00113
00114 std::vector<std::pair<const PSimHit*, int > >::const_iterator simHitIter = input.begin();
00115 std::vector<std::pair<const PSimHit*, int > >::const_iterator simHitIterEnd = input.end();
00116 if(theFlatDistribution->fire()>inefficiency) {
00117 for (;simHitIter != simHitIterEnd; ++simHitIter) {
00118
00119 if ( std::fabs( ((*simHitIter).first)->tof() - cosmicShift - det->surface().toGlobal(((*simHitIter).first)->localPosition()).mag()/30.) < tofCut && ((*simHitIter).first)->energyLoss()>0) {
00120 localFirstChannel = numStrips;
00121 localLastChannel = 0;
00122
00123 theSiHitDigitizer->processHit(((*simHitIter).first),*det,bfield,langle, locAmpl, localFirstChannel, localLastChannel);
00124
00125
00126
00127
00128 if(APVSaturationFromHIP&&!zeroSuppression){
00129 int pdg_id = ((*simHitIter).first)->particleType();
00130 particle = pdt->particle(pdg_id);
00131 if(particle != NULL){
00132 float charge = particle->charge();
00133 bool isHadron = particle->isHadron();
00134 if(charge!=0 && isHadron){
00135 if(theFlatDistribution->fire()<APVSaturationProb){
00136 int FirstAPV = localFirstChannel/128;
00137 int LastAPV = localLastChannel/128;
00138 std::cout << "-------------------HIP--------------" << std::endl;
00139 std::cout << "Killing APVs " << FirstAPV << " - " <<LastAPV << " " << detID <<std::endl;
00140 for(int strip = FirstAPV*128; strip < LastAPV*128 +128; ++strip) badChannels[strip] = true;
00141
00142
00143
00144 }
00145 }
00146 }
00147 }
00148
00149
00150 theSiPileUpSignals->add(locAmpl, localFirstChannel, localLastChannel, ((*simHitIter).first), (*simHitIter).second);
00151
00152
00153 for (size_t iChannel=localFirstChannel; iChannel<localLastChannel; iChannel++) {
00154 if(locAmpl[iChannel]>0.) {
00155
00156
00157 detAmpl[iChannel]+=locAmpl[iChannel];
00158 }
00159 }
00160 if(firstChannelWithSignal>localFirstChannel) firstChannelWithSignal=localFirstChannel;
00161 if(lastChannelWithSignal<localLastChannel) lastChannelWithSignal=localLastChannel;
00162 }
00163 }
00164 }
00165
00166
00167 for(int strip =0; strip < numStrips; ++strip) if(badChannels[strip]) detAmpl[strip] =0;
00168
00169 const SiPileUpSignals::HitToDigisMapType& theLink(theSiPileUpSignals->dumpLink());
00170 const SiPileUpSignals::HitCounterToDigisMapType& theCounterLink(theSiPileUpSignals->dumpCounterLink());
00171
00172 if(zeroSuppression){
00173 if(noise){
00174 int RefStrip = 0;
00175 while(RefStrip<numStrips&&badChannels[RefStrip]){
00176 RefStrip++;
00177 }
00178 if(RefStrip<numStrips){
00179 float noiseRMS = noiseHandle->getNoise(RefStrip,detNoiseRange);
00180 float gainValue = gainHandle->getStripGain(RefStrip, detGainRange);
00181 theSiNoiseAdder->addNoise(detAmpl,firstChannelWithSignal,lastChannelWithSignal,numStrips,noiseRMS*theElectronPerADC/gainValue);
00182 }
00183 }
00184 digis.clear();
00185 theSiZeroSuppress->suppress(theSiDigitalConverter->convert(detAmpl, gainHandle, detID), digis, detID,noiseHandle,thresholdHandle);
00186 push_link(digis, theLink, theCounterLink, detAmpl,detID);
00187 outdigi.data = digis;
00188 }
00189
00190
00191 if(!zeroSuppression){
00192
00193
00194
00195
00196
00197
00198
00199
00200
00201
00202
00203
00204
00205
00206
00207
00208
00209
00210 if(BaselineShift){
00211 theSiNoiseAdder->addBaselineShift(detAmpl, badChannels);
00212 }
00213
00214
00215
00216 if(noise){
00217 std::vector<float> noiseRMSv;
00218 noiseRMSv.clear();
00219 noiseRMSv.insert(noiseRMSv.begin(),numStrips,0.);
00220
00221 if(SingleStripNoise){
00222 for(int strip=0; strip< numStrips; ++strip){
00223 if(!badChannels[strip]) noiseRMSv[strip] = (noiseHandle->getNoise(strip,detNoiseRange))* theElectronPerADC;
00224 }
00225
00226 } else {
00227 int RefStrip = 0;
00228 while(RefStrip<numStrips&&badChannels[RefStrip]){
00229 RefStrip++;
00230 }
00231 if(RefStrip<numStrips){
00232 float noiseRMS = noiseHandle->getNoise(RefStrip,detNoiseRange) *theElectronPerADC;
00233 for(int strip=0; strip< numStrips; ++strip){
00234 if(!badChannels[strip]) noiseRMSv[strip] = noiseRMS;
00235 }
00236 }
00237 }
00238
00239 theSiNoiseAdder->addNoiseVR(detAmpl, noiseRMSv);
00240 }
00241
00242
00243
00244 if(CommonModeNoise){
00245 float cmnRMS = 0.;
00246 DetId detId(detID);
00247 uint32_t SubDet = detId.subdetId();
00248 if(SubDet==3){
00249 cmnRMS = cmnRMStib;
00250 }else if(SubDet==4){
00251 cmnRMS = cmnRMStid;
00252 }else if(SubDet==5){
00253 cmnRMS = cmnRMStob;
00254 }else if(SubDet==6){
00255 cmnRMS = cmnRMStec;
00256 }
00257 cmnRMS *= theElectronPerADC;
00258 theSiNoiseAdder->addCMNoise(detAmpl, cmnRMS, badChannels);
00259 }
00260
00261
00262
00263
00264
00265 std::vector<float> vPeds;
00266 vPeds.clear();
00267 vPeds.insert(vPeds.begin(),numStrips,0.);
00268
00269 if(RealPedestals){
00270 for(int strip=0; strip< numStrips; ++strip){
00271 if(!badChannels[strip]) vPeds[strip] = (pedestalHandle->getPed(strip,detPedestalRange)+pedOffset)* theElectronPerADC;
00272 }
00273 } else {
00274 for(int strip=0; strip< numStrips; ++strip){
00275 if(!badChannels[strip]) vPeds[strip] = pedOffset* theElectronPerADC;
00276 }
00277 }
00278
00279 theSiNoiseAdder->addPedestals(detAmpl, vPeds);
00280
00281
00282
00283
00284
00285
00286 rawdigis.clear();
00287 rawdigis = theSiDigitalConverter->convertRaw(detAmpl, gainHandle, detID);
00288 push_link_raw(rawdigis, theLink, theCounterLink, detAmpl,detID);
00289 outrawdigi.data = rawdigis;
00290
00291
00292 }
00293 }
00294
00295 void SiStripDigitizerAlgorithm::push_link(const DigitalVecType &digis,
00296 const HitToDigisMapType& htd,
00297 const HitCounterToDigisMapType& hctd,
00298 const std::vector<double>& afterNoise,
00299 unsigned int detID) {
00300 link_coll.clear();
00301 for ( DigitalVecType::const_iterator i=digis.begin(); i!=digis.end(); i++) {
00302
00303
00304 HitToDigisMapType::const_iterator mi(htd.find(i->strip()));
00305 if (mi == htd.end()) continue;
00306 HitCounterToDigisMapType::const_iterator cmi(hctd.find(i->strip()));
00307 std::map<const PSimHit *, Amplitude> totalAmplitudePerSimHit;
00308 for (std::vector < std::pair < const PSimHit*, Amplitude > >::const_iterator simul =
00309 (*mi).second.begin() ; simul != (*mi).second.end(); simul ++){
00310 totalAmplitudePerSimHit[(*simul).first] += (*simul).second;
00311 }
00312
00313
00314 double totalAmplitude1 = afterNoise[(*mi).first];
00315
00316
00317 int sim_counter=0;
00318 for (std::map<const PSimHit *, Amplitude>::const_iterator iter = totalAmplitudePerSimHit.begin();
00319 iter != totalAmplitudePerSimHit.end(); iter++){
00320 float threshold = 0.;
00321 float fraction = (*iter).second/totalAmplitude1;
00322 if ( fraction >= threshold) {
00323
00324 if(fraction > 1.) fraction = 1.;
00325 for (std::vector < std::pair < const PSimHit*, int > >::const_iterator
00326 simcount = (*cmi).second.begin() ; simcount != (*cmi).second.end(); simcount ++){
00327 if((*iter).first == (*simcount).first) sim_counter = (*simcount).second;
00328 }
00329 link_coll.push_back(StripDigiSimLink( (*mi).first,
00330 ((*iter).first)->trackId(),
00331 sim_counter,
00332 ((*iter).first)->eventId(),
00333 fraction));
00334 }
00335 }
00336 }
00337 }
00338
00339 void SiStripDigitizerAlgorithm::push_link_raw(const DigitalRawVecType &digis,
00340 const HitToDigisMapType& htd,
00341 const HitCounterToDigisMapType& hctd,
00342 const std::vector<double>& afterNoise,
00343 unsigned int detID) {
00344 link_coll.clear();
00345 int nstrip = -1;
00346 for ( DigitalRawVecType::const_iterator i=digis.begin(); i!=digis.end(); i++) {
00347 nstrip++;
00348
00349
00350 HitToDigisMapType::const_iterator mi(htd.find(nstrip));
00351 HitCounterToDigisMapType::const_iterator cmi(hctd.find(nstrip));
00352 if (mi == htd.end()) continue;
00353 std::map<const PSimHit *, Amplitude> totalAmplitudePerSimHit;
00354 for (std::vector < std::pair < const PSimHit*, Amplitude > >::const_iterator simul =
00355 (*mi).second.begin() ; simul != (*mi).second.end(); simul ++){
00356 totalAmplitudePerSimHit[(*simul).first] += (*simul).second;
00357 }
00358
00359
00360 double totalAmplitude1 = afterNoise[(*mi).first];
00361
00362
00363 int sim_counter_raw=0;
00364 for (std::map<const PSimHit *, Amplitude>::const_iterator iter = totalAmplitudePerSimHit.begin();
00365 iter != totalAmplitudePerSimHit.end(); iter++){
00366 float threshold = 0.;
00367 float fraction = (*iter).second/totalAmplitude1;
00368 if (fraction >= threshold) {
00369
00370 if(fraction >1.) fraction = 1.;
00371
00372 for (std::vector < std::pair < const PSimHit*, int > >::const_iterator
00373 simcount = (*cmi).second.begin() ; simcount != (*cmi).second.end(); simcount ++){
00374 if((*iter).first == (*simcount).first) sim_counter_raw = (*simcount).second;
00375 }
00376 link_coll.push_back(StripDigiSimLink( (*mi).first,
00377 ((*iter).first)->trackId(),
00378 sim_counter_raw,
00379 ((*iter).first)->eventId(),
00380 fraction));
00381 }
00382 }
00383 }
00384 }