CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SiStripNoiseNormalizedWithApvGainBuilder.cc
Go to the documentation of this file.
10 #include <iostream>
11 #include <fstream>
12 #include <vector>
13 #include <algorithm>
14 
16  fp_(iConfig.getUntrackedParameter<edm::FileInPath>("file",edm::FileInPath("CalibTracker/SiStripCommon/data/SiStripDetInfo.dat"))),
17  printdebug_(iConfig.getUntrackedParameter<uint32_t>("printDebug",1)),
18  pset_(iConfig),
19  electronsPerADC_(0.),
20  minimumPosValue_(0.),
21  stripLengthMode_(true),
22  printDebug_(0)
23 {}
24 
26 {
27  // Read the gain from the given tag
28  edm::ESHandle<SiStripApvGain> inputApvGain;
29  iSetup.get<SiStripApvGainRcd>().get( inputApvGain );
30  std::vector<uint32_t> inputDetIds;
31  inputApvGain->getDetIds(inputDetIds);
32 
33  // Prepare the new object
35 
37 
38 
39 
40  stripLengthMode_ = pset_.getParameter<bool>("StripLengthMode");
41 
42  //parameters for random noise generation. not used if Strip length mode is chosen
43  std::map<int, std::vector<double> > meanNoise;
44  fillParameters(meanNoise, "MeanNoise");
45  std::map<int, std::vector<double> > sigmaNoise;
46  fillParameters(sigmaNoise, "SigmaNoise");
47  minimumPosValue_ = pset_.getParameter<double>("MinPositiveNoise");
48 
49  //parameters for strip length proportional noise generation. not used if random mode is chosen
50  std::map<int, std::vector<double> > noiseStripLengthLinearSlope;
51  fillParameters(noiseStripLengthLinearSlope, "NoiseStripLengthSlope");
52  std::map<int, std::vector<double> > noiseStripLengthLinearQuote;
53  fillParameters(noiseStripLengthLinearQuote, "NoiseStripLengthQuote");
54  electronsPerADC_ = pset_.getParameter<double>("electronPerAdc");
55 
56  printDebug_ = pset_.getUntrackedParameter<uint32_t>("printDebug", 5);
57 
58  unsigned int count = 0;
59  const std::map<uint32_t, SiStripDetInfoFileReader::DetInfo > DetInfos = reader.getAllData();
60  for(std::map<uint32_t, SiStripDetInfoFileReader::DetInfo >::const_iterator it = DetInfos.begin(); it != DetInfos.end(); it++) {
61 
62  // Find if this DetId is in the input tag and if so how many are the Apvs for which it contains information
63  SiStripApvGain::Range inputRange(inputApvGain->getRange(it->first));
64 
65  //Generate Noises for det detid
66  SiStripNoises::InputVector theSiStripVector;
67  float noise = 0.;
68  uint32_t detId = it->first;
69  std::pair<int, int> sl = subDetAndLayer(detId);
70  unsigned short nApvs = it->second.nApvs;
71 
72  if(stripLengthMode_) {
73  // Use strip length
74  double linearSlope = noiseStripLengthLinearSlope[sl.first][sl.second];
75  double linearQuote = noiseStripLengthLinearQuote[sl.first][sl.second];
76  double stripLength = it->second.stripLength;
77  for( unsigned short j=0; j<nApvs; ++j ) {
78 
79  double gain = inputApvGain->getApvGain(j, inputRange);
80 
81  for( unsigned short stripId = 0; stripId < 128; ++stripId ) {
82  noise = ( ( linearSlope*stripLength + linearQuote) / electronsPerADC_ ) * gain;
83  if( count<printDebug_ ) printLog(detId, stripId+128*j, noise);
84  obj->setData(noise, theSiStripVector);
85  }
86  }
87  }
88  else {
89  // Use random generator
90  double meanN = meanNoise[sl.first][sl.second];
91  double sigmaN = sigmaNoise[sl.first][sl.second];
92  for( unsigned short j=0; j<nApvs; ++j ) {
93 
94  double gain = inputApvGain->getApvGain(j, inputRange);
95 
96  for( unsigned short stripId = 0; stripId < 128; ++stripId ) {
97  noise = ( CLHEP::RandGauss::shoot(meanN, sigmaN) ) * gain;
98  if( noise<=minimumPosValue_ ) noise = minimumPosValue_;
99  if( count<printDebug_ ) printLog(detId, stripId+128*j, noise);
100  obj->setData(noise, theSiStripVector);
101  }
102  }
103  }
104  ++count;
105 
106  if ( ! obj->put(it->first,theSiStripVector) ) {
107  edm::LogError("SiStripNoisesFakeESSource::produce ")<<" detid already exists"<<std::endl;
108  }
109  }
110 
111  //End now write data in DB
113 
114  if( mydbservice.isAvailable() ){
115  if( mydbservice->isNewTagRequest("SiStripNoisesRcd") ){
116  mydbservice->createNewIOV<SiStripNoises>(obj,mydbservice->beginOfTime(),mydbservice->endOfTime(),"SiStripNoisesRcd");
117  }
118  else {
119  mydbservice->appendSinceTime<SiStripNoises>(obj,mydbservice->currentTime(),"SiStripNoisesRcd");
120  }
121  }
122  else {
123  edm::LogError("SiStripNoiseNormalizedWithApvGainBuilder")<<"Service is unavailable"<<std::endl;
124  }
125 }
126 
127 std::pair<int, int> SiStripNoiseNormalizedWithApvGainBuilder::subDetAndLayer( const uint32_t detId ) const
128 {
129  int layerId = 0;
130 
131  StripSubdetector subid(detId);
132  int subId = subid.subdetId();
133 
134  if( subId == int(StripSubdetector::TIB)) {
135  TIBDetId theTIBDetId(detId);
136  layerId = theTIBDetId.layer() - 1;
137  }
138  else if(subId == int(StripSubdetector::TOB)) {
139  TOBDetId theTOBDetId(detId);
140  layerId = theTOBDetId.layer() - 1;
141  }
142  else if(subId == int(StripSubdetector::TID)) {
143  TIDDetId theTIDDetId(detId);
144  layerId = theTIDDetId.ring() - 1;
145  }
146  if(subId == int(StripSubdetector::TEC)) {
147  TECDetId theTECDetId = TECDetId(detId);
148  layerId = theTECDetId.ring() - 1;
149  }
150  return std::make_pair(subId, layerId);
151 }
152 
153 void SiStripNoiseNormalizedWithApvGainBuilder::fillParameters(std::map<int, std::vector<double> > & mapToFill, const std::string & parameterName) const
154 {
155  int layersTIB = 4;
156  int ringsTID = 3;
157  int layersTOB = 6;
158  int ringsTEC = 7;
159 
160  fillSubDetParameter( mapToFill, pset_.getParameter<std::vector<double> >(parameterName+"TIB"), int(StripSubdetector::TIB), layersTIB );
161  fillSubDetParameter( mapToFill, pset_.getParameter<std::vector<double> >(parameterName+"TID"), int(StripSubdetector::TID), ringsTID );
162  fillSubDetParameter( mapToFill, pset_.getParameter<std::vector<double> >(parameterName+"TOB"), int(StripSubdetector::TOB), layersTOB );
163  fillSubDetParameter( mapToFill, pset_.getParameter<std::vector<double> >(parameterName+"TEC"), int(StripSubdetector::TEC), ringsTEC );
164 }
165 
166 void SiStripNoiseNormalizedWithApvGainBuilder::fillSubDetParameter(std::map<int, std::vector<double> > & mapToFill, const std::vector<double> & v, const int subDet, const unsigned short layers) const
167 {
168  if( v.size() == layers ) {
169  mapToFill.insert(std::make_pair( subDet, v ));
170  }
171  else if( v.size() == 1 ) {
172  std::vector<double> parV(layers, v[0]);
173  mapToFill.insert(std::make_pair( subDet, parV ));
174  }
175  else {
176  throw cms::Exception("Configuration") << "ERROR: number of parameters for subDet " << subDet << " are " << v.size() << ". They must be either 1 or " << layers << std::endl;
177  }
178 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
inputRange
Get input source.
unsigned int layer() const
layer id
Definition: TOBDetId.h:39
std::vector< uint16_t > InputVector
Definition: SiStripNoises.h:44
SiStripNoiseNormalizedWithApvGainBuilder(const edm::ParameterSet &iConfig)
void appendSinceTime(T *payloadObj, cond::Time_t sinceTime, const std::string &recordName, bool withlogging=false)
bool isNewTagRequest(const std::string &recordName)
unsigned int ring() const
ring id
Definition: TIDDetId.h:55
bool isAvailable() const
Definition: Service.h:47
std::pair< ContainerIterator, ContainerIterator > Range
int j
Definition: DBlmapReader.cc:9
void fillParameters(std::map< int, std::vector< double > > &mapToFill, const std::string &parameterName) const
Fills the parameters read from cfg and matching the name in the given map.
void printLog(const uint32_t detId, const unsigned short strip, const double &noise) const
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:39
void createNewIOV(T *firstPayloadObj, cond::Time_t firstSinceTime, cond::Time_t firstTillTime, const std::string &recordName, bool withlogging=false)
bool put(const uint32_t &detID, const InputVector &input)
void fillSubDetParameter(std::map< int, std::vector< double > > &mapToFill, const std::vector< double > &v, const int subDet, const unsigned short layers) const
virtual void analyze(const edm::Event &, const edm::EventSetup &)
const T & get() const
Definition: EventSetup.h:55
unsigned int layer() const
layer id
Definition: TIBDetId.h:41
unsigned int ring() const
ring id
Definition: TECDetId.h:71
std::string fullPath() const
Definition: FileInPath.cc:171
std::pair< int, int > subDetAndLayer(const uint32_t detit) const
Given the map and the detid it returns the corresponding layer/ring.
mathSSE::Vec4< T > v
void setData(float noise_, InputVector &vped)