Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 11:59:51

0001 #include "CalibTracker/SiStripESProducers/plugins/fake/SiStripNoiseNormalizedWithApvGainBuilder.h"
0002 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
0003 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0004 #include "Geometry/TrackerGeometryBuilder/interface/StripGeomDetUnit.h"
0005 #include "FWCore/Framework/interface/EventSetup.h"
0006 #include "SiStripFakeAPVParameters.h"
0007 #include <iostream>
0008 #include <fstream>
0009 #include <vector>
0010 #include <algorithm>
0011 
0012 SiStripNoiseNormalizedWithApvGainBuilder::SiStripNoiseNormalizedWithApvGainBuilder(const edm::ParameterSet& iConfig)
0013     : printdebug_(iConfig.getUntrackedParameter<uint32_t>("printDebug", 1)),
0014       pset_(iConfig),
0015       electronsPerADC_(0.),
0016       minimumPosValue_(0.),
0017       stripLengthMode_(true),
0018       printDebug_(0) {
0019   tTopoToken_ = esConsumes<TrackerTopology, TrackerTopologyRcd>();
0020   tGeomToken_ = esConsumes<TrackerGeometry, TrackerDigiGeometryRecord>();
0021   inputApvGainToken_ = esConsumes<SiStripApvGain, SiStripApvGainRcd>();
0022 }
0023 
0024 void SiStripNoiseNormalizedWithApvGainBuilder::analyze(const edm::Event& evt, const edm::EventSetup& iSetup) {
0025   const auto& tTopo = iSetup.getData(tTopoToken_);
0026   const auto& tGeom = iSetup.getData(tGeomToken_);
0027   // Read the gain from the given tag
0028   const auto& inputApvGain = iSetup.getData(inputApvGainToken_);
0029   std::vector<uint32_t> inputDetIds;
0030   inputApvGain.getDetIds(inputDetIds);
0031 
0032   // Prepare the new object
0033   SiStripNoises obj;
0034 
0035   stripLengthMode_ = pset_.getParameter<bool>("StripLengthMode");
0036 
0037   //parameters for random noise generation. not used if Strip length mode is chosen
0038   SiStripFakeAPVParameters meanNoise{pset_, "MeanNoise"};
0039   SiStripFakeAPVParameters sigmaNoise{pset_, "SigmaNoise"};
0040   minimumPosValue_ = pset_.getParameter<double>("MinPositiveNoise");
0041 
0042   //parameters for strip length proportional noise generation. not used if random mode is chosen
0043   SiStripFakeAPVParameters noiseStripLengthLinearSlope{pset_, "NoiseStripLengthSlope"};
0044   SiStripFakeAPVParameters noiseStripLengthLinearQuote{pset_, "NoiseStripLengthQuote"};
0045   electronsPerADC_ = pset_.getParameter<double>("electronPerAdc");
0046 
0047   printDebug_ = pset_.getUntrackedParameter<uint32_t>("printDebug", 5);
0048 
0049   unsigned int count = 0;
0050   for (const auto det : tGeom.detUnits()) {
0051     const StripGeomDetUnit* stripDet = dynamic_cast<const StripGeomDetUnit*>(det);
0052     if (stripDet != nullptr) {
0053       const DetId detId = stripDet->geographicalId();
0054       // Find if this DetId is in the input tag and if so how many are the Apvs for which it contains information
0055       SiStripApvGain::Range inputRange(inputApvGain.getRange(detId));
0056 
0057       //Generate Noises for det detid
0058       SiStripNoises::InputVector theSiStripVector;
0059       float noise = 0.;
0060       SiStripFakeAPVParameters::index sl = SiStripFakeAPVParameters::getIndex(&tTopo, detId);
0061       unsigned short nApvs = stripDet->specificTopology().nstrips() / 128;
0062 
0063       if (stripLengthMode_) {
0064         // Use strip length
0065         double linearSlope = noiseStripLengthLinearSlope.get(sl);
0066         double linearQuote = noiseStripLengthLinearQuote.get(sl);
0067         double stripLength = stripDet->specificTopology().stripLength();
0068         for (unsigned short j = 0; j < nApvs; ++j) {
0069           double gain = inputApvGain.getApvGain(j, inputRange);
0070 
0071           for (unsigned short stripId = 0; stripId < 128; ++stripId) {
0072             noise = ((linearSlope * stripLength + linearQuote) / electronsPerADC_) * gain;
0073             if (count < printDebug_)
0074               printLog(detId, stripId + 128 * j, noise);
0075             obj.setData(noise, theSiStripVector);
0076           }
0077         }
0078       } else {
0079         // Use random generator
0080         double meanN = meanNoise.get(sl);
0081         double sigmaN = sigmaNoise.get(sl);
0082         for (unsigned short j = 0; j < nApvs; ++j) {
0083           double gain = inputApvGain.getApvGain(j, inputRange);
0084 
0085           for (unsigned short stripId = 0; stripId < 128; ++stripId) {
0086             noise = (CLHEP::RandGauss::shoot(meanN, sigmaN)) * gain;
0087             if (noise <= minimumPosValue_)
0088               noise = minimumPosValue_;
0089             if (count < printDebug_)
0090               printLog(detId, stripId + 128 * j, noise);
0091             obj.setData(noise, theSiStripVector);
0092           }
0093         }
0094       }
0095       ++count;
0096 
0097       if (!obj.put(detId, theSiStripVector)) {
0098         edm::LogError("SiStripNoisesFakeESSource::produce ") << " detid already exists" << std::endl;
0099       }
0100     }
0101   }
0102 
0103   //End now write data in DB
0104   edm::Service<cond::service::PoolDBOutputService> mydbservice;
0105 
0106   if (mydbservice.isAvailable()) {
0107     if (mydbservice->isNewTagRequest("SiStripNoisesRcd")) {
0108       mydbservice->createOneIOV<SiStripNoises>(obj, mydbservice->beginOfTime(), "SiStripNoisesRcd");
0109     } else {
0110       mydbservice->appendOneIOV<SiStripNoises>(obj, mydbservice->currentTime(), "SiStripNoisesRcd");
0111     }
0112   } else {
0113     edm::LogError("SiStripNoiseNormalizedWithApvGainBuilder") << "Service is unavailable" << std::endl;
0114   }
0115 }