CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_14/src/SimTracker/SiStripDigitizer/src/SiStripDigitizerAlgorithm.cc

Go to the documentation of this file.
00001 // File: SiStripDigitizerAlgorithm.cc
00002 // Description:  Steering class for digitization.
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   //delete particle;
00067   //delete pdt;
00068 }
00069 
00070 //  Run the algorithm for a given module
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   //stroing the bad stip of the the module. the module is not removed but just signal put to 0
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   // local amplitude of detector channels (from processed PSimHit)
00104 //  locAmpl.clear();
00105   detAmpl.clear();
00106 //  locAmpl.insert(locAmpl.begin(),numStrips,0.);
00107   // total amplitude of detector channels
00108   detAmpl.insert(detAmpl.begin(),numStrips,0.);
00109 
00110   firstChannelWithSignal = numStrips;
00111   lastChannelWithSignal  = 0;
00112 
00113   // First: loop on the SimHits
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       locAmpl.clear();
00119       locAmpl.insert(locAmpl.begin(),numStrips,0.);
00120       // check TOF
00121       if ( std::fabs( ((*simHitIter).first)->tof() - cosmicShift - det->surface().toGlobal(((*simHitIter).first)->localPosition()).mag()/30.) < tofCut && ((*simHitIter).first)->energyLoss()>0) {
00122         localFirstChannel = numStrips;
00123         localLastChannel  = 0;
00124         // process the hit
00125         theSiHitDigitizer->processHit(((*simHitIter).first),*det,bfield,langle, locAmpl, localFirstChannel, localLastChannel);
00126           
00127                   //APV Killer to simulate HIP effect
00128                   //------------------------------------------------------
00129                   
00130                   if(APVSaturationFromHIP&&!zeroSuppression){
00131                     int pdg_id = ((*simHitIter).first)->particleType();
00132                         particle = pdt->particle(pdg_id);
00133                         if(particle != NULL){
00134                                 float charge = particle->charge();
00135                                 bool isHadron = particle->isHadron();
00136                             if(charge!=0 && isHadron){
00137                                         if(theFlatDistribution->fire()<APVSaturationProb){
00138                                                 int FirstAPV = localFirstChannel/128;
00139                                                 int LastAPV = localLastChannel/128;
00140                                                 std::cout << "-------------------HIP--------------" << std::endl;
00141                                                 std::cout << "Killing APVs " << FirstAPV << " - " <<LastAPV << " " << detID <<std::endl;
00142                                                 for(int strip = FirstAPV*128; strip < LastAPV*128 +128; ++strip) badChannels[strip] = true; //doing like that I remove the signal information only after the 
00143                                                                                                                                                                                                                                         //stip that got the HIP but it remains the signal of the previous
00144                                                                                                                                                                                                                                         //one. I'll make a further loop to remove all signal                                                                                                                                                                                                            
00145                                                 
00146                                         }
00147                                 }
00148                         }
00149               }             
00150                 
00151                 
00152                 theSiPileUpSignals->add(locAmpl, localFirstChannel, localLastChannel, ((*simHitIter).first), (*simHitIter).second);
00153     
00154                 // sum signal on strips
00155         for (size_t iChannel=localFirstChannel; iChannel<localLastChannel; iChannel++) {
00156           if(locAmpl[iChannel]>0.) {
00157                     //if(!badChannels[iChannel]) detAmpl[iChannel]+=locAmpl[iChannel];
00158                         //locAmpl[iChannel]=0;
00159                         detAmpl[iChannel]+=locAmpl[iChannel];
00160            }
00161         }
00162         if(firstChannelWithSignal>localFirstChannel) firstChannelWithSignal=localFirstChannel;
00163         if(lastChannelWithSignal<localLastChannel) lastChannelWithSignal=localLastChannel;
00164       }
00165         }
00166   }
00167   
00168   //removing signal from the dead (and HIP effected) strips
00169   for(int strip =0; strip < numStrips; ++strip) if(badChannels[strip]) detAmpl[strip] =0;
00170   
00171   const SiPileUpSignals::HitToDigisMapType& theLink(theSiPileUpSignals->dumpLink());  
00172   const SiPileUpSignals::HitCounterToDigisMapType& theCounterLink(theSiPileUpSignals->dumpCounterLink());  
00173   
00174   if(zeroSuppression){
00175     if(noise){
00176           int RefStrip = int(numStrips/2.);
00177           while(RefStrip<numStrips&&badChannels[RefStrip]){ //if the refstrip is bad, I move up to when I don't find it
00178                 RefStrip++;
00179           }
00180           if(RefStrip<numStrips){
00181                 float noiseRMS = noiseHandle->getNoise(RefStrip,detNoiseRange);
00182                 float gainValue = gainHandle->getStripGain(RefStrip, detGainRange);
00183                 theSiNoiseAdder->addNoise(detAmpl,firstChannelWithSignal,lastChannelWithSignal,numStrips,noiseRMS*theElectronPerADC/gainValue);
00184           }
00185         }
00186     digis.clear();
00187     theSiZeroSuppress->suppress(theSiDigitalConverter->convert(detAmpl, gainHandle, detID), digis, detID,noiseHandle,thresholdHandle);
00188     push_link(digis, theLink, theCounterLink, detAmpl,detID);
00189     outdigi.data = digis;
00190   }
00191   
00192   
00193   if(!zeroSuppression){
00194     //if(noise){
00195       // the constant pedestal offset is needed because
00196       //   negative adc counts are not allowed in case
00197       //   Pedestal and CMN subtraction is performed.
00198       //   The pedestal value read from the conditions
00199       //   is pedValue and after the pedestal subtraction
00200       //   the baseline is zero. The Common Mode Noise
00201       //   is not subtracted from the negative adc counts
00202       //   channels. Adding pedOffset the baseline is set
00203       //   to pedOffset after pedestal subtraction and CMN
00204       //   is subtracted to all the channels since none of
00205       //   them has negative adc value. The pedOffset is
00206       //   treated as a constant component in the CMN
00207       //   estimation and subtracted as CMN.
00208       
00209          
00210                 //calculating the charge deposited on each APV and subtracting the shift
00211                 //------------------------------------------------------
00212                 if(BaselineShift){
00213                    theSiNoiseAdder->addBaselineShift(detAmpl, badChannels);
00214                 }
00215                 
00216                 //Adding the strip noise
00217                 //------------------------------------------------------                                                 
00218                 if(noise){
00219                     std::vector<float> noiseRMSv;
00220                         noiseRMSv.clear();
00221                     noiseRMSv.insert(noiseRMSv.begin(),numStrips,0.);
00222                         
00223                     if(SingleStripNoise){
00224                             for(int strip=0; strip< numStrips; ++strip){
00225                                         if(!badChannels[strip]) noiseRMSv[strip] = (noiseHandle->getNoise(strip,detNoiseRange))* theElectronPerADC;
00226                                 }
00227                         
00228                 } else {
00229                             int RefStrip = 0; //int(numStrips/2.);
00230                             while(RefStrip<numStrips&&badChannels[RefStrip]){ //if the refstrip is bad, I move up to when I don't find it
00231                                         RefStrip++;
00232                                 }
00233                                 if(RefStrip<numStrips){
00234                                         float noiseRMS = noiseHandle->getNoise(RefStrip,detNoiseRange) *theElectronPerADC;
00235                                         for(int strip=0; strip< numStrips; ++strip){
00236                                         if(!badChannels[strip]) noiseRMSv[strip] = noiseRMS;
00237                                         }
00238                                 }
00239                         }
00240                         
00241                         theSiNoiseAdder->addNoiseVR(detAmpl, noiseRMSv);
00242                 }                       
00243                 
00244                 //adding the CMN
00245                 //------------------------------------------------------
00246         if(CommonModeNoise){
00247                   float cmnRMS = 0.;
00248                   DetId  detId(detID);
00249                   uint32_t SubDet = detId.subdetId();
00250                   if(SubDet==3){
00251                     cmnRMS = cmnRMStib;
00252                   }else if(SubDet==4){
00253                     cmnRMS = cmnRMStid;
00254                   }else if(SubDet==5){
00255                     cmnRMS = cmnRMStob;
00256                   }else if(SubDet==6){
00257                     cmnRMS = cmnRMStec;
00258                   }
00259                   cmnRMS *= theElectronPerADC;
00260           theSiNoiseAdder->addCMNoise(detAmpl, cmnRMS, badChannels);
00261                 }
00262                 
00263                         
00264                 //Adding the pedestals
00265                 //------------------------------------------------------
00266                 
00267                 std::vector<float> vPeds;
00268                 vPeds.clear();
00269                 vPeds.insert(vPeds.begin(),numStrips,0.);
00270                 
00271                 if(RealPedestals){
00272                     for(int strip=0; strip< numStrips; ++strip){
00273                            if(!badChannels[strip]) vPeds[strip] = (pedestalHandle->getPed(strip,detPedestalRange)+pedOffset)* theElectronPerADC;
00274                     }
00275         } else {
00276                     for(int strip=0; strip< numStrips; ++strip){
00277                           if(!badChannels[strip]) vPeds[strip] = pedOffset* theElectronPerADC;
00278                         }
00279                 }
00280                 
00281                 theSiNoiseAdder->addPedestals(detAmpl, vPeds);  
00282                 
00283                  
00284         //if(!RealPedestals&&!CommonModeNoise&&!noise&&!BaselineShift&&!APVSaturationFromHIP){
00285     //  edm::LogWarning("SiStripDigitizer")<<"You are running the digitizer without Noise generation and without applying Zero Suppression. ARE YOU SURE???";
00286     //}else{                                                     
00287     
00288         rawdigis.clear();
00289     rawdigis = theSiDigitalConverter->convertRaw(detAmpl, gainHandle, detID);
00290     push_link_raw(rawdigis, theLink, theCounterLink, detAmpl,detID);
00291     outrawdigi.data = rawdigis;
00292         
00293         //}
00294   }
00295 }
00296 
00297 void SiStripDigitizerAlgorithm::push_link(const DigitalVecType &digis,
00298                                           const HitToDigisMapType& htd,
00299                                           const HitCounterToDigisMapType& hctd,
00300                                           const std::vector<double>& afterNoise,
00301                                           unsigned int detID) {
00302   link_coll.clear();  
00303   for ( DigitalVecType::const_iterator i=digis.begin(); i!=digis.end(); i++) {
00304     // Instead of checking the validity of the links against the digis,
00305     //  let's loop over digis and push the corresponding link
00306     HitToDigisMapType::const_iterator mi(htd.find(i->strip()));  
00307     if (mi == htd.end()) continue;
00308     HitCounterToDigisMapType::const_iterator cmi(hctd.find(i->strip()));  
00309     std::map<const PSimHit *, Amplitude> totalAmplitudePerSimHit;
00310     for (std::vector < std::pair < const PSimHit*, Amplitude > >::const_iterator simul = 
00311            (*mi).second.begin() ; simul != (*mi).second.end(); simul ++){
00312       totalAmplitudePerSimHit[(*simul).first] += (*simul).second;
00313     }
00314     
00315     //--- include the noise as well
00316     double totalAmplitude1 = afterNoise[(*mi).first];
00317     
00318     //--- digisimlink
00319     int sim_counter=0; 
00320     for (std::map<const PSimHit *, Amplitude>::const_iterator iter = totalAmplitudePerSimHit.begin(); 
00321          iter != totalAmplitudePerSimHit.end(); iter++){
00322       float threshold = 0.;
00323       float fraction = (*iter).second/totalAmplitude1;
00324       if ( fraction >= threshold) {
00325         // Noise fluctuation could make fraction>1. Unphysical, set it by hand = 1.
00326         if(fraction > 1.) fraction = 1.;
00327         for (std::vector < std::pair < const PSimHit*, int > >::const_iterator 
00328                simcount = (*cmi).second.begin() ; simcount != (*cmi).second.end(); simcount ++){
00329           if((*iter).first == (*simcount).first) sim_counter = (*simcount).second;
00330         }
00331         link_coll.push_back(StripDigiSimLink( (*mi).first, //channel
00332                                               ((*iter).first)->trackId(), //simhit trackId
00333                                               sim_counter, //simhit counter
00334                                               ((*iter).first)->eventId(), //simhit eventId
00335                                               fraction)); //fraction
00336       }
00337     }
00338   }
00339 }
00340 
00341 void SiStripDigitizerAlgorithm::push_link_raw(const DigitalRawVecType &digis,
00342                                               const HitToDigisMapType& htd,
00343                                               const HitCounterToDigisMapType& hctd,
00344                                               const std::vector<double>& afterNoise,
00345                                               unsigned int detID) {
00346   link_coll.clear();  
00347   int nstrip = -1;
00348   for ( DigitalRawVecType::const_iterator i=digis.begin(); i!=digis.end(); i++) {
00349     nstrip++;
00350     // Instead of checking the validity of the links against the digis,
00351     //  let's loop over digis and push the corresponding link
00352     HitToDigisMapType::const_iterator mi(htd.find(nstrip));  
00353     HitCounterToDigisMapType::const_iterator cmi(hctd.find(nstrip));  
00354     if (mi == htd.end()) continue;
00355     std::map<const PSimHit *, Amplitude> totalAmplitudePerSimHit;
00356     for (std::vector < std::pair < const PSimHit*, Amplitude > >::const_iterator simul = 
00357            (*mi).second.begin() ; simul != (*mi).second.end(); simul ++){
00358       totalAmplitudePerSimHit[(*simul).first] += (*simul).second;
00359     }
00360     
00361     //--- include the noise as well
00362     double totalAmplitude1 = afterNoise[(*mi).first];
00363     
00364     //--- digisimlink
00365     int sim_counter_raw=0;
00366     for (std::map<const PSimHit *, Amplitude>::const_iterator iter = totalAmplitudePerSimHit.begin(); 
00367          iter != totalAmplitudePerSimHit.end(); iter++){
00368       float threshold = 0.;
00369       float fraction = (*iter).second/totalAmplitude1;
00370       if (fraction >= threshold) {
00371         //Noise fluctuation could make fraction>1. Unphysical, set it by hand.
00372         if(fraction >1.) fraction = 1.;
00373         //add counter information
00374         for (std::vector < std::pair < const PSimHit*, int > >::const_iterator 
00375                simcount = (*cmi).second.begin() ; simcount != (*cmi).second.end(); simcount ++){
00376           if((*iter).first == (*simcount).first) sim_counter_raw = (*simcount).second;
00377         }
00378         link_coll.push_back(StripDigiSimLink( (*mi).first, //channel
00379                                               ((*iter).first)->trackId(), //simhit trackId
00380                                               sim_counter_raw, //simhit counter
00381                                               ((*iter).first)->eventId(), //simhit eventId
00382                                               fraction)); //fraction
00383       }
00384     }
00385   }
00386 }