CMS 3D CMS Logo

SiStripApvGainBuilderFromTag.cc
Go to the documentation of this file.
7 #include <iostream>
8 #include <fstream>
9 #include <vector>
10 #include <algorithm>
11 
13 
15  printdebug_(iConfig.getUntrackedParameter<uint32_t>("printDebug",1)),
16  pset_(iConfig)
17 {}
18 
20 {
21  //Retrieve tracker topology from geometry
23  iSetup.get<IdealGeometryRecord>().get(tTopoHandle);
24  const TrackerTopology* const tTopo = tTopoHandle.product();
25 
26  // unsigned int run=evt.id().run();
27 
28  std::string genMode = pset_.getParameter<std::string>("genMode");
29  bool applyTuning = pset_.getParameter<bool>("applyTuning");
30 
31  double meanGain_=pset_.getParameter<double>("MeanGain");
32  double sigmaGain_=pset_.getParameter<double>("SigmaGain");
33  double minimumPosValue_=pset_.getParameter<double>("MinPositiveGain");
34 
35  uint32_t printdebug_ = pset_.getUntrackedParameter<uint32_t>("printDebug", 5);
36 
37  //parameters for layer/disk level correction; not used if applyTuning=false
38  SiStripFakeAPVParameters correct{pset_, "correct"};
39 
40  // Read the gain from the given tag
41  edm::ESHandle<SiStripApvGain> inputApvGain;
42  iSetup.get<SiStripApvGainRcd>().get( inputApvGain );
43  std::vector<uint32_t> inputDetIds;
44  inputApvGain->getDetIds(inputDetIds);
45 
46  // Prepare the new object
48 
50  uint32_t count = 0;
51  const std::map<uint32_t, SiStripDetInfoFileReader::DetInfo >& DetInfos = reader->getAllData();
52  for(std::map<uint32_t, SiStripDetInfoFileReader::DetInfo >::const_iterator it = DetInfos.begin(); it != DetInfos.end(); it++) {
53 
54  // Find if this DetId is in the input tag and if so how many are the Apvs for which it contains information
56  size_t inputRangeSize = 0;
57  if( find( inputDetIds.begin(), inputDetIds.end(), it->first ) != inputDetIds.end() ) {
58  inputRange = inputApvGain->getRange(it->first);
59  inputRangeSize = distance(inputRange.first, inputRange.second);
60  }
61 
62  std::vector<float> theSiStripVector;
63  for(unsigned short j=0; j<it->second.nApvs; j++){
64 
65  double gainValue = meanGain_;
66 
67  if( j < inputRangeSize ) {
68  gainValue = inputApvGain->getApvGain(j, inputRange);
69  // cout << "Gain = " << gainValue <<" from input tag for DetId = " << it->first << " and apv = " << j << endl;
70  }
71  // else {
72  // cout << "No gain in input tag for DetId = " << it->first << " and apv = " << j << " using value from cfg = " << gainValue << endl;
73  // }
74 
75 
76  // corrections at layer/disk level:
77  uint32_t detId = it->first;
79  //unsigned short nApvs = it->second.nApvs;
80  if (applyTuning) {
81  double correction = correct.get(sl);
82  gainValue *= correction;
83  }
84 
85  // smearing:
86  if (genMode == "gaussian") {
87  gainValue = CLHEP::RandGauss::shoot(gainValue, sigmaGain_);
88  if(gainValue<=minimumPosValue_) gainValue=minimumPosValue_;
89  }
90  else if( genMode != "default" ) {
91  LogDebug("SiStripApvGain") << "ERROR: wrong genMode specifier : " << genMode << ", please select one of \"default\" or \"gaussian\"" << std::endl;
92  exit(1);
93  }
94 
95  if (count<printdebug_) {
96  edm::LogInfo("SiStripApvGainGeneratorFromTag") << "detid: " << it->first << " Apv: " << j << " gain: " << gainValue << std::endl;
97  }
98  theSiStripVector.push_back(gainValue);
99  }
100  count++;
101  SiStripApvGain::Range range(theSiStripVector.begin(),theSiStripVector.end());
102  if ( ! obj->put(it->first,range) )
103  edm::LogError("SiStripApvGainGeneratorFromTag")<<" detid already exists"<<std::endl;
104  }
105 
106  //End now write data in DB
108 
109  if( mydbservice.isAvailable() ){
110  if( mydbservice->isNewTagRequest("SiStripApvGainRcd2") ){
111  mydbservice->createNewIOV<SiStripApvGain>(obj,mydbservice->beginOfTime(),mydbservice->endOfTime(),"SiStripApvGainRcd2");
112  }
113  else {
114  mydbservice->appendSinceTime<SiStripApvGain>(obj,mydbservice->currentTime(),"SiStripApvGainRcd2");
115  }
116  }
117  else {
118  edm::LogError("SiStripApvGainBuilderFromTag")<<"Service is unavailable"<<std::endl;
119  }
120 }
#define LogDebug(id)
T getParameter(std::string const &) const
SiStripApvGainBuilderFromTag(const edm::ParameterSet &iConfig)
T getUntrackedParameter(std::string const &, T const &) const
static float getApvGain(uint16_t apv, const Range &range)
static index getIndex(const TrackerTopology *tTopo, DetId id)
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
void getDetIds(std::vector< uint32_t > &DetIds_) const
bool put(const uint32_t &detID, Range input)
const std::map< uint32_t, DetInfo > & getAllData() const
void appendSinceTime(T *payloadObj, cond::Time_t sinceTime, const std::string &recordName, bool withlogging=false)
bool isNewTagRequest(const std::string &recordName)
bool isAvailable() const
Definition: Service.h:40
std::pair< ContainerIterator, ContainerIterator > Range
void createNewIOV(T *firstPayloadObj, cond::Time_t firstSinceTime, cond::Time_t firstTillTime, const std::string &recordName, bool withlogging=false)
inputRange
Get input source.
T get() const
Definition: EventSetup.h:71
const Range getRange(const uint32_t detID) const
void analyze(const edm::Event &, const edm::EventSetup &) override
T const * product() const
Definition: ESHandle.h:86