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