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/SiStripApvGainBuilderFromTag.h"
0002 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
0003 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0004 #include "Geometry/TrackerGeometryBuilder/interface/StripGeomDetUnit.h"
0005 #include <iostream>
0006 #include <fstream>
0007 #include <vector>
0008 #include <algorithm>
0009 
0010 #include "SiStripFakeAPVParameters.h"
0011 
0012 SiStripApvGainBuilderFromTag::SiStripApvGainBuilderFromTag(const edm::ParameterSet& iConfig)
0013     : printdebug_(iConfig.getUntrackedParameter<uint32_t>("printDebug", 1)), pset_(iConfig) {
0014   tTopoToken_ = esConsumes<TrackerTopology, TrackerTopologyRcd>();
0015   tGeomToken_ = esConsumes<TrackerGeometry, TrackerDigiGeometryRecord>();
0016   inputApvGainToken_ = esConsumes<SiStripApvGain, SiStripApvGainRcd>();
0017 }
0018 
0019 void SiStripApvGainBuilderFromTag::analyze(const edm::Event& evt, const edm::EventSetup& iSetup) {
0020   const auto& tTopo = iSetup.getData(tTopoToken_);
0021   const auto& tGeom = iSetup.getData(tGeomToken_);
0022 
0023   //   unsigned int run=evt.id().run();
0024 
0025   std::string genMode = pset_.getParameter<std::string>("genMode");
0026   bool applyTuning = pset_.getParameter<bool>("applyTuning");
0027 
0028   double meanGain_ = pset_.getParameter<double>("MeanGain");
0029   double sigmaGain_ = pset_.getParameter<double>("SigmaGain");
0030   double minimumPosValue_ = pset_.getParameter<double>("MinPositiveGain");
0031 
0032   uint32_t printdebug_ = pset_.getUntrackedParameter<uint32_t>("printDebug", 5);
0033 
0034   //parameters for layer/disk level correction; not used if applyTuning=false
0035   SiStripFakeAPVParameters correct{pset_, "correct"};
0036 
0037   // Read the gain from the given tag
0038   const auto& inputApvGain = iSetup.getData(inputApvGainToken_);
0039   std::vector<uint32_t> inputDetIds;
0040   inputApvGain.getDetIds(inputDetIds);
0041 
0042   // Prepare the new object
0043   SiStripApvGain obj;
0044 
0045   uint32_t count = 0;
0046   for (const auto det : tGeom.detUnits()) {
0047     const StripGeomDetUnit* stripDet = dynamic_cast<const StripGeomDetUnit*>(det);
0048     if (stripDet != nullptr) {
0049       const DetId detid = stripDet->geographicalId();
0050       // Find if this DetId is in the input tag and if so how many are the Apvs for which it contains information
0051       SiStripApvGain::Range inputRange;
0052       size_t inputRangeSize = 0;
0053       if (find(inputDetIds.begin(), inputDetIds.end(), detid) != inputDetIds.end()) {
0054         inputRange = inputApvGain.getRange(detid);
0055         inputRangeSize = distance(inputRange.first, inputRange.second);
0056       }
0057 
0058       std::vector<float> theSiStripVector;
0059       for (unsigned short j = 0; j < (stripDet->specificTopology().nstrips() / 128); j++) {
0060         double gainValue = meanGain_;
0061 
0062         if (j < inputRangeSize) {
0063           gainValue = inputApvGain.getApvGain(j, inputRange);
0064           // cout << "Gain = " << gainValue <<" from input tag for DetId = " << detid.rawId() << " and apv = " << j << endl;
0065         }
0066         // else {
0067         //   cout << "No gain in input tag for DetId = " << detid << " and apv = " << j << " using value from cfg = " << gainValue << endl;
0068         // }
0069 
0070         // corrections at layer/disk level:
0071         SiStripFakeAPVParameters::index sl = SiStripFakeAPVParameters::getIndex(&tTopo, detid);
0072         //unsigned short nApvs = (stripDet->specificTopology().nstrips()/128);
0073         if (applyTuning) {
0074           double correction = correct.get(sl);
0075           gainValue *= correction;
0076         }
0077 
0078         // smearing:
0079         if (genMode == "gaussian") {
0080           gainValue = CLHEP::RandGauss::shoot(gainValue, sigmaGain_);
0081           if (gainValue <= minimumPosValue_)
0082             gainValue = minimumPosValue_;
0083         } else if (genMode != "default") {
0084           LogDebug("SiStripApvGain") << "ERROR: wrong genMode specifier : " << genMode
0085                                      << ", please select one of \"default\" or \"gaussian\"" << std::endl;
0086           exit(1);
0087         }
0088 
0089         if (count < printdebug_) {
0090           edm::LogInfo("SiStripApvGainGeneratorFromTag")
0091               << "detid: " << detid.rawId() << " Apv: " << j << " gain: " << gainValue << std::endl;
0092         }
0093         theSiStripVector.push_back(gainValue);
0094       }
0095       count++;
0096       SiStripApvGain::Range range(theSiStripVector.begin(), theSiStripVector.end());
0097       if (!obj.put(detid, range))
0098         edm::LogError("SiStripApvGainGeneratorFromTag") << " detid already exists" << std::endl;
0099     }
0100   }
0101 
0102   //End now write data in DB
0103   edm::Service<cond::service::PoolDBOutputService> mydbservice;
0104 
0105   if (mydbservice.isAvailable()) {
0106     if (mydbservice->isNewTagRequest("SiStripApvGainRcd2")) {
0107       mydbservice->createOneIOV<SiStripApvGain>(obj, mydbservice->beginOfTime(), "SiStripApvGainRcd2");
0108     } else {
0109       mydbservice->appendOneIOV<SiStripApvGain>(obj, mydbservice->currentTime(), "SiStripApvGainRcd2");
0110     }
0111   } else {
0112     edm::LogError("SiStripApvGainBuilderFromTag") << "Service is unavailable" << std::endl;
0113   }
0114 }