CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_14/src/CalibTracker/SiStripESProducers/plugins/fake/SiStripNoiseNormalizedWithApvGainBuilder.cc

Go to the documentation of this file.
00001 #include "CalibTracker/SiStripESProducers/plugins/fake/SiStripNoiseNormalizedWithApvGainBuilder.h"
00002 #include "CalibTracker/SiStripCommon/interface/SiStripDetInfoFileReader.h"
00003 #include "CondFormats/DataRecord/interface/SiStripCondDataRecords.h"
00004 #include "DataFormats/SiStripDetId/interface/StripSubdetector.h"
00005 #include "DataFormats/SiStripDetId/interface/TIBDetId.h"
00006 #include "DataFormats/SiStripDetId/interface/TIDDetId.h"
00007 #include "DataFormats/SiStripDetId/interface/TOBDetId.h"
00008 #include "DataFormats/SiStripDetId/interface/TECDetId.h"
00009 #include "FWCore/Framework/interface/EventSetup.h"
00010 #include <iostream>
00011 #include <fstream>
00012 #include <vector>
00013 #include <algorithm>
00014 
00015 SiStripNoiseNormalizedWithApvGainBuilder::SiStripNoiseNormalizedWithApvGainBuilder( const edm::ParameterSet& iConfig ):
00016   fp_(iConfig.getUntrackedParameter<edm::FileInPath>("file",edm::FileInPath("CalibTracker/SiStripCommon/data/SiStripDetInfo.dat"))),
00017   printdebug_(iConfig.getUntrackedParameter<uint32_t>("printDebug",1)),
00018   pset_(iConfig),
00019   electronsPerADC_(0.),
00020   minimumPosValue_(0.),
00021   stripLengthMode_(true),
00022   printDebug_(0)
00023 {}
00024 
00025 void SiStripNoiseNormalizedWithApvGainBuilder::analyze(const edm::Event& evt, const edm::EventSetup& iSetup)
00026 {
00027   // Read the gain from the given tag
00028   edm::ESHandle<SiStripApvGain> inputApvGain;
00029   iSetup.get<SiStripApvGainRcd>().get( inputApvGain );
00030   std::vector<uint32_t> inputDetIds;
00031   inputApvGain->getDetIds(inputDetIds);
00032 
00033   // Prepare the new object
00034   SiStripNoises* obj = new SiStripNoises();
00035 
00036   SiStripDetInfoFileReader reader(fp_.fullPath());
00037 
00038 
00039 
00040   stripLengthMode_ = pset_.getParameter<bool>("StripLengthMode");
00041   
00042   //parameters for random noise generation. not used if Strip length mode is chosen
00043   std::map<int, std::vector<double> > meanNoise;
00044   fillParameters(meanNoise, "MeanNoise");
00045   std::map<int, std::vector<double> > sigmaNoise;
00046   fillParameters(sigmaNoise, "SigmaNoise");
00047   minimumPosValue_ = pset_.getParameter<double>("MinPositiveNoise");
00048 
00049   //parameters for strip length proportional noise generation. not used if random mode is chosen
00050   std::map<int, std::vector<double> > noiseStripLengthLinearSlope;
00051   fillParameters(noiseStripLengthLinearSlope, "NoiseStripLengthSlope");
00052   std::map<int, std::vector<double> > noiseStripLengthLinearQuote;
00053   fillParameters(noiseStripLengthLinearQuote, "NoiseStripLengthQuote");
00054   electronsPerADC_ = pset_.getParameter<double>("electronPerAdc");        
00055 
00056   printDebug_ = pset_.getUntrackedParameter<uint32_t>("printDebug", 5);
00057 
00058   unsigned int count = 0;
00059   const std::map<uint32_t, SiStripDetInfoFileReader::DetInfo > DetInfos = reader.getAllData();
00060   for(std::map<uint32_t, SiStripDetInfoFileReader::DetInfo >::const_iterator it = DetInfos.begin(); it != DetInfos.end(); it++) {
00061 
00062     // Find if this DetId is in the input tag and if so how many are the Apvs for which it contains information
00063     SiStripApvGain::Range inputRange(inputApvGain->getRange(it->first));
00064 
00065     //Generate Noises for det detid
00066     SiStripNoises::InputVector theSiStripVector;
00067     float noise = 0.;
00068     uint32_t detId = it->first;
00069     std::pair<int, int> sl = subDetAndLayer(detId);
00070     unsigned short nApvs = it->second.nApvs;
00071 
00072     if(stripLengthMode_) {
00073       // Use strip length
00074       double linearSlope = noiseStripLengthLinearSlope[sl.first][sl.second];
00075       double linearQuote = noiseStripLengthLinearQuote[sl.first][sl.second];
00076       double stripLength = it->second.stripLength;
00077       for( unsigned short j=0; j<nApvs; ++j ) {
00078 
00079         double gain = inputApvGain->getApvGain(j, inputRange);
00080 
00081         for( unsigned short stripId = 0; stripId < 128; ++stripId ) {
00082           noise = ( ( linearSlope*stripLength + linearQuote) / electronsPerADC_ ) * gain;
00083           if( count<printDebug_ ) printLog(detId, stripId+128*j, noise);
00084           obj->setData(noise, theSiStripVector);
00085         }
00086       }
00087     }
00088     else {
00089       // Use random generator
00090       double meanN = meanNoise[sl.first][sl.second];
00091       double sigmaN = sigmaNoise[sl.first][sl.second];
00092       for( unsigned short j=0; j<nApvs; ++j ) {
00093 
00094         double gain = inputApvGain->getApvGain(j, inputRange);
00095 
00096         for( unsigned short stripId = 0; stripId < 128; ++stripId ) {
00097           noise = ( CLHEP::RandGauss::shoot(meanN, sigmaN) ) * gain;
00098           if( noise<=minimumPosValue_ ) noise = minimumPosValue_;
00099           if( count<printDebug_ ) printLog(detId, stripId+128*j, noise);
00100           obj->setData(noise, theSiStripVector);
00101         }
00102       }
00103     }
00104     ++count;
00105     
00106     if ( ! obj->put(it->first,theSiStripVector) ) {
00107       edm::LogError("SiStripNoisesFakeESSource::produce ")<<" detid already exists"<<std::endl;
00108     }
00109   }
00110   
00111   //End now write data in DB
00112   edm::Service<cond::service::PoolDBOutputService> mydbservice;
00113   
00114   if( mydbservice.isAvailable() ){
00115     if( mydbservice->isNewTagRequest("SiStripNoisesRcd") ){
00116       mydbservice->createNewIOV<SiStripNoises>(obj,mydbservice->beginOfTime(),mydbservice->endOfTime(),"SiStripNoisesRcd");
00117     }
00118     else {
00119       mydbservice->appendSinceTime<SiStripNoises>(obj,mydbservice->currentTime(),"SiStripNoisesRcd");
00120     }
00121   }
00122   else {
00123     edm::LogError("SiStripNoiseNormalizedWithApvGainBuilder")<<"Service is unavailable"<<std::endl;
00124   }
00125 }
00126 
00127 std::pair<int, int> SiStripNoiseNormalizedWithApvGainBuilder::subDetAndLayer( const uint32_t detId ) const
00128 {
00129   int layerId = 0;
00130 
00131   StripSubdetector subid(detId);
00132   int subId = subid.subdetId();
00133 
00134   if( subId == int(StripSubdetector::TIB)) {
00135     TIBDetId theTIBDetId(detId);
00136     layerId = theTIBDetId.layer() - 1;
00137   }
00138   else if(subId == int(StripSubdetector::TOB)) {
00139     TOBDetId theTOBDetId(detId);
00140     layerId = theTOBDetId.layer() - 1;
00141   }
00142   else if(subId == int(StripSubdetector::TID)) {
00143     TIDDetId theTIDDetId(detId);
00144     layerId = theTIDDetId.ring() - 1;
00145   }
00146   if(subId == int(StripSubdetector::TEC)) {
00147     TECDetId theTECDetId = TECDetId(detId); 
00148     layerId = theTECDetId.ring() - 1;
00149   }
00150   return std::make_pair(subId, layerId);
00151 }
00152 
00153 void SiStripNoiseNormalizedWithApvGainBuilder::fillParameters(std::map<int, std::vector<double> > & mapToFill, const std::string & parameterName) const
00154 {
00155   int layersTIB = 4;
00156   int ringsTID = 3;
00157   int layersTOB = 6;
00158   int ringsTEC = 7;
00159 
00160   fillSubDetParameter( mapToFill, pset_.getParameter<std::vector<double> >(parameterName+"TIB"), int(StripSubdetector::TIB), layersTIB );
00161   fillSubDetParameter( mapToFill, pset_.getParameter<std::vector<double> >(parameterName+"TID"), int(StripSubdetector::TID), ringsTID );
00162   fillSubDetParameter( mapToFill, pset_.getParameter<std::vector<double> >(parameterName+"TOB"), int(StripSubdetector::TOB), layersTOB );
00163   fillSubDetParameter( mapToFill, pset_.getParameter<std::vector<double> >(parameterName+"TEC"), int(StripSubdetector::TEC), ringsTEC );
00164 }
00165 
00166 void SiStripNoiseNormalizedWithApvGainBuilder::fillSubDetParameter(std::map<int, std::vector<double> > & mapToFill, const std::vector<double> & v, const int subDet, const unsigned short layers) const
00167 {
00168   if( v.size() == layers ) {
00169     mapToFill.insert(std::make_pair( subDet, v ));
00170   }
00171   else if( v.size() == 1 ) {
00172     std::vector<double> parV(layers, v[0]);
00173     mapToFill.insert(std::make_pair( subDet, parV ));
00174   }
00175   else {
00176     throw cms::Exception("Configuration") << "ERROR: number of parameters for subDet " << subDet << " are " << v.size() << ". They must be either 1 or " << layers << std::endl;
00177   }
00178 }