CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
SiStripNoiseNormalizedWithApvGainBuilder Class Reference

#include <SiStripNoiseNormalizedWithApvGainBuilder.h>

Inheritance diagram for SiStripNoiseNormalizedWithApvGainBuilder:
edm::EDAnalyzer edm::EDConsumerBase

Public Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &)
 
 SiStripNoiseNormalizedWithApvGainBuilder (const edm::ParameterSet &iConfig)
 
 ~SiStripNoiseNormalizedWithApvGainBuilder ()
 
- Public Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzer ()
 
ModuleDescription const & moduleDescription () const
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 
- Public Member Functions inherited from edm::EDConsumerBase
 EDConsumerBase ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesDependentUpon (const std::string &iProcessName, std::vector< const char * > &oModuleLabels) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Private Member Functions

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. More...
 
void fillSubDetParameter (std::map< int, std::vector< double > > &mapToFill, const std::vector< double > &v, const int subDet, const unsigned short layers) const
 
void printLog (const uint32_t detId, const unsigned short strip, const double &noise) const
 
std::pair< int, int > subDetAndLayer (const uint32_t detit, const TrackerTopology *tTopo) const
 Given the map and the detid it returns the corresponding layer/ring. More...
 

Private Attributes

double electronsPerADC_
 
edm::FileInPath fp_
 
double minimumPosValue_
 
bool printdebug_
 
uint32_t printDebug_
 
edm::ParameterSet pset_
 
bool stripLengthMode_
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Produces a noise tag using the same settings as the service used in the DummyDBWriter, but it receives a SiStripApvGain tag from the EventSetup and uses the gain values (per apv) to rescale the noise (per strip).

Definition at line 29 of file SiStripNoiseNormalizedWithApvGainBuilder.h.

Constructor & Destructor Documentation

SiStripNoiseNormalizedWithApvGainBuilder::SiStripNoiseNormalizedWithApvGainBuilder ( const edm::ParameterSet iConfig)
explicit
SiStripNoiseNormalizedWithApvGainBuilder::~SiStripNoiseNormalizedWithApvGainBuilder ( )
inline

Definition at line 35 of file SiStripNoiseNormalizedWithApvGainBuilder.h.

35 {};

Member Function Documentation

void SiStripNoiseNormalizedWithApvGainBuilder::analyze ( const edm::Event evt,
const edm::EventSetup iSetup 
)
virtual

Implements edm::EDAnalyzer.

Definition at line 23 of file SiStripNoiseNormalizedWithApvGainBuilder.cc.

References cond::service::PoolDBOutputService::appendSinceTime(), cond::service::PoolDBOutputService::beginOfTime(), prof2calltree::count, cond::service::PoolDBOutputService::createNewIOV(), cond::service::PoolDBOutputService::currentTime(), electronsPerADC_, cond::service::PoolDBOutputService::endOfTime(), fillParameters(), fp_, edm::FileInPath::fullPath(), edm::EventSetup::get(), edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), estimatePileup::inputRange, edm::Service< T >::isAvailable(), cond::service::PoolDBOutputService::isNewTagRequest(), j, minimumPosValue_, getGTfromDQMFile::obj, printDebug_, printLog(), edm::ESHandle< class >::product(), pset_, SiStripNoises::put(), matplotRender::reader, SiStripNoises::setData(), stripLengthMode_, and subDetAndLayer().

24 {
25  //Retrieve tracker topology from geometry
27  iSetup.get<IdealGeometryRecord>().get(tTopoHandle);
28  const TrackerTopology* const tTopo = tTopoHandle.product();
29 
30  // Read the gain from the given tag
31  edm::ESHandle<SiStripApvGain> inputApvGain;
32  iSetup.get<SiStripApvGainRcd>().get( inputApvGain );
33  std::vector<uint32_t> inputDetIds;
34  inputApvGain->getDetIds(inputDetIds);
35 
36  // Prepare the new object
38 
40 
41 
42 
43  stripLengthMode_ = pset_.getParameter<bool>("StripLengthMode");
44 
45  //parameters for random noise generation. not used if Strip length mode is chosen
46  std::map<int, std::vector<double> > meanNoise;
47  fillParameters(meanNoise, "MeanNoise");
48  std::map<int, std::vector<double> > sigmaNoise;
49  fillParameters(sigmaNoise, "SigmaNoise");
50  minimumPosValue_ = pset_.getParameter<double>("MinPositiveNoise");
51 
52  //parameters for strip length proportional noise generation. not used if random mode is chosen
53  std::map<int, std::vector<double> > noiseStripLengthLinearSlope;
54  fillParameters(noiseStripLengthLinearSlope, "NoiseStripLengthSlope");
55  std::map<int, std::vector<double> > noiseStripLengthLinearQuote;
56  fillParameters(noiseStripLengthLinearQuote, "NoiseStripLengthQuote");
57  electronsPerADC_ = pset_.getParameter<double>("electronPerAdc");
58 
59  printDebug_ = pset_.getUntrackedParameter<uint32_t>("printDebug", 5);
60 
61  unsigned int count = 0;
62  const std::map<uint32_t, SiStripDetInfoFileReader::DetInfo > DetInfos = reader.getAllData();
63  for(std::map<uint32_t, SiStripDetInfoFileReader::DetInfo >::const_iterator it = DetInfos.begin(); it != DetInfos.end(); it++) {
64 
65  // Find if this DetId is in the input tag and if so how many are the Apvs for which it contains information
66  SiStripApvGain::Range inputRange(inputApvGain->getRange(it->first));
67 
68  //Generate Noises for det detid
69  SiStripNoises::InputVector theSiStripVector;
70  float noise = 0.;
71  uint32_t detId = it->first;
72  std::pair<int, int> sl = subDetAndLayer(detId, tTopo);
73  unsigned short nApvs = it->second.nApvs;
74 
75  if(stripLengthMode_) {
76  // Use strip length
77  double linearSlope = noiseStripLengthLinearSlope[sl.first][sl.second];
78  double linearQuote = noiseStripLengthLinearQuote[sl.first][sl.second];
79  double stripLength = it->second.stripLength;
80  for( unsigned short j=0; j<nApvs; ++j ) {
81 
82  double gain = inputApvGain->getApvGain(j, inputRange);
83 
84  for( unsigned short stripId = 0; stripId < 128; ++stripId ) {
85  noise = ( ( linearSlope*stripLength + linearQuote) / electronsPerADC_ ) * gain;
86  if( count<printDebug_ ) printLog(detId, stripId+128*j, noise);
87  obj->setData(noise, theSiStripVector);
88  }
89  }
90  }
91  else {
92  // Use random generator
93  double meanN = meanNoise[sl.first][sl.second];
94  double sigmaN = sigmaNoise[sl.first][sl.second];
95  for( unsigned short j=0; j<nApvs; ++j ) {
96 
97  double gain = inputApvGain->getApvGain(j, inputRange);
98 
99  for( unsigned short stripId = 0; stripId < 128; ++stripId ) {
100  noise = ( CLHEP::RandGauss::shoot(meanN, sigmaN) ) * gain;
101  if( noise<=minimumPosValue_ ) noise = minimumPosValue_;
102  if( count<printDebug_ ) printLog(detId, stripId+128*j, noise);
103  obj->setData(noise, theSiStripVector);
104  }
105  }
106  }
107  ++count;
108 
109  if ( ! obj->put(it->first,theSiStripVector) ) {
110  edm::LogError("SiStripNoisesFakeESSource::produce ")<<" detid already exists"<<std::endl;
111  }
112  }
113 
114  //End now write data in DB
116 
117  if( mydbservice.isAvailable() ){
118  if( mydbservice->isNewTagRequest("SiStripNoisesRcd") ){
119  mydbservice->createNewIOV<SiStripNoises>(obj,mydbservice->beginOfTime(),mydbservice->endOfTime(),"SiStripNoisesRcd");
120  }
121  else {
122  mydbservice->appendSinceTime<SiStripNoises>(obj,mydbservice->currentTime(),"SiStripNoisesRcd");
123  }
124  }
125  else {
126  edm::LogError("SiStripNoiseNormalizedWithApvGainBuilder")<<"Service is unavailable"<<std::endl;
127  }
128 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
inputRange
Get input source.
std::pair< int, int > subDetAndLayer(const uint32_t detit, const TrackerTopology *tTopo) const
Given the map and the detid it returns the corresponding layer/ring.
std::vector< uint16_t > InputVector
Definition: SiStripNoises.h:51
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:46
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
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)
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:62
std::string fullPath() const
Definition: FileInPath.cc:165
void setData(float noise_, InputVector &vped)
void SiStripNoiseNormalizedWithApvGainBuilder::fillParameters ( std::map< int, std::vector< double > > &  mapToFill,
const std::string &  parameterName 
) const
private

Fills the parameters read from cfg and matching the name in the given map.

Definition at line 152 of file SiStripNoiseNormalizedWithApvGainBuilder.cc.

References fillSubDetParameter(), edm::ParameterSet::getParameter(), pset_, StripSubdetector::TEC, StripSubdetector::TIB, StripSubdetector::TID, and StripSubdetector::TOB.

Referenced by analyze().

153 {
154  int layersTIB = 4;
155  int ringsTID = 3;
156  int layersTOB = 6;
157  int ringsTEC = 7;
158 
159  fillSubDetParameter( mapToFill, pset_.getParameter<std::vector<double> >(parameterName+"TIB"), int(StripSubdetector::TIB), layersTIB );
160  fillSubDetParameter( mapToFill, pset_.getParameter<std::vector<double> >(parameterName+"TID"), int(StripSubdetector::TID), ringsTID );
161  fillSubDetParameter( mapToFill, pset_.getParameter<std::vector<double> >(parameterName+"TOB"), int(StripSubdetector::TOB), layersTOB );
162  fillSubDetParameter( mapToFill, pset_.getParameter<std::vector<double> >(parameterName+"TEC"), int(StripSubdetector::TEC), ringsTEC );
163 }
T getParameter(std::string const &) const
void fillSubDetParameter(std::map< int, std::vector< double > > &mapToFill, const std::vector< double > &v, const int subDet, const unsigned short layers) const
void SiStripNoiseNormalizedWithApvGainBuilder::fillSubDetParameter ( std::map< int, std::vector< double > > &  mapToFill,
const std::vector< double > &  v,
const int  subDet,
const unsigned short  layers 
) const
private

Fills the map with the paramters for the given subdetector.
Each vector "v" holds the parameters for the layers/rings, if the vector has only one parameter all the layers/rings get that parameter.
The only other possibility is that the number of parameters equals the number of layers, otherwise an exception of type "Configuration" will be thrown.

Definition at line 165 of file SiStripNoiseNormalizedWithApvGainBuilder.cc.

References edm::hlt::Exception, and LayerTriplets::layers().

Referenced by fillParameters().

166 {
167  if( v.size() == layers ) {
168  mapToFill.insert(std::make_pair( subDet, v ));
169  }
170  else if( v.size() == 1 ) {
171  std::vector<double> parV(layers, v[0]);
172  mapToFill.insert(std::make_pair( subDet, parV ));
173  }
174  else {
175  throw cms::Exception("Configuration") << "ERROR: number of parameters for subDet " << subDet << " are " << v.size() << ". They must be either 1 or " << layers << std::endl;
176  }
177 }
std::vector< LayerSetAndLayers > layers(const SeedingLayerSetsHits &sets)
Definition: LayerTriplets.cc:4
void SiStripNoiseNormalizedWithApvGainBuilder::printLog ( const uint32_t  detId,
const unsigned short  strip,
const double &  noise 
) const
inlineprivate

Definition at line 53 of file SiStripNoiseNormalizedWithApvGainBuilder.h.

Referenced by analyze().

54  {
55  edm::LogInfo("SiStripNoisesDummyCalculator") << "detid: " << detId << " strip: " << strip << " noise: " << noise << " \t" << std::endl;
56  }
std::pair< int, int > SiStripNoiseNormalizedWithApvGainBuilder::subDetAndLayer ( const uint32_t  detit,
const TrackerTopology tTopo 
) const
private

Given the map and the detid it returns the corresponding layer/ring.

Definition at line 130 of file SiStripNoiseNormalizedWithApvGainBuilder.cc.

References DetId::subdetId(), StripSubdetector::TEC, TrackerTopology::tecRing(), StripSubdetector::TIB, TrackerTopology::tibLayer(), StripSubdetector::TID, TrackerTopology::tidRing(), StripSubdetector::TOB, and TrackerTopology::tobLayer().

Referenced by analyze().

131 {
132  int layerId = 0;
133 
134  StripSubdetector subid(detId);
135  int subId = subid.subdetId();
136 
137  if( subId == int(StripSubdetector::TIB)) {
138  layerId = tTopo->tibLayer(detId) - 1;
139  }
140  else if(subId == int(StripSubdetector::TOB)) {
141  layerId = tTopo->tobLayer(detId) - 1;
142  }
143  else if(subId == int(StripSubdetector::TID)) {
144  layerId = tTopo->tidRing(detId) - 1;
145  }
146  if(subId == int(StripSubdetector::TEC)) {
147  layerId = tTopo->tecRing(detId) - 1;
148  }
149  return std::make_pair(subId, layerId);
150 }
unsigned int tibLayer(const DetId &id) const
unsigned int tidRing(const DetId &id) const
unsigned int tecRing(const DetId &id) const
ring id
unsigned int tobLayer(const DetId &id) const

Member Data Documentation

double SiStripNoiseNormalizedWithApvGainBuilder::electronsPerADC_
private

Definition at line 62 of file SiStripNoiseNormalizedWithApvGainBuilder.h.

Referenced by analyze().

edm::FileInPath SiStripNoiseNormalizedWithApvGainBuilder::fp_
private

Definition at line 58 of file SiStripNoiseNormalizedWithApvGainBuilder.h.

Referenced by analyze().

double SiStripNoiseNormalizedWithApvGainBuilder::minimumPosValue_
private

Definition at line 63 of file SiStripNoiseNormalizedWithApvGainBuilder.h.

Referenced by analyze().

bool SiStripNoiseNormalizedWithApvGainBuilder::printdebug_
private

Definition at line 59 of file SiStripNoiseNormalizedWithApvGainBuilder.h.

uint32_t SiStripNoiseNormalizedWithApvGainBuilder::printDebug_
private

Definition at line 65 of file SiStripNoiseNormalizedWithApvGainBuilder.h.

Referenced by analyze().

edm::ParameterSet SiStripNoiseNormalizedWithApvGainBuilder::pset_
private

Definition at line 60 of file SiStripNoiseNormalizedWithApvGainBuilder.h.

Referenced by analyze(), and fillParameters().

bool SiStripNoiseNormalizedWithApvGainBuilder::stripLengthMode_
private

Definition at line 64 of file SiStripNoiseNormalizedWithApvGainBuilder.h.

Referenced by analyze().