CMS 3D CMS Logo

Public Member Functions | Private Member Functions | Private Attributes

SiStripNoiseNormalizedWithApvGainBuilder Class Reference

#include <SiStripNoiseNormalizedWithApvGainBuilder.h>

Inheritance diagram for SiStripNoiseNormalizedWithApvGainBuilder:
edm::EDAnalyzer

List of all members.

Public Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &)
 SiStripNoiseNormalizedWithApvGainBuilder (const edm::ParameterSet &iConfig)
 ~SiStripNoiseNormalizedWithApvGainBuilder ()

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.
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
 Given the map and the detid it returns the corresponding layer/ring.

Private Attributes

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

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 28 of file SiStripNoiseNormalizedWithApvGainBuilder.h.


Constructor & Destructor Documentation

SiStripNoiseNormalizedWithApvGainBuilder::SiStripNoiseNormalizedWithApvGainBuilder ( const edm::ParameterSet iConfig) [explicit]

Definition at line 15 of file SiStripNoiseNormalizedWithApvGainBuilder.cc.

                                                                                                                  :
  fp_(iConfig.getUntrackedParameter<edm::FileInPath>("file",edm::FileInPath("CalibTracker/SiStripCommon/data/SiStripDetInfo.dat"))),
  printdebug_(iConfig.getUntrackedParameter<uint32_t>("printDebug",1)),
  pset_(iConfig),
  electronsPerADC_(0.),
  minimumPosValue_(0.),
  stripLengthMode_(true),
  printDebug_(0)
{}
SiStripNoiseNormalizedWithApvGainBuilder::~SiStripNoiseNormalizedWithApvGainBuilder ( ) [inline]

Definition at line 34 of file SiStripNoiseNormalizedWithApvGainBuilder.h.

{};

Member Function Documentation

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

Implements edm::EDAnalyzer.

Definition at line 25 of file SiStripNoiseNormalizedWithApvGainBuilder.cc.

References prof2calltree::count, electronsPerADC_, fillParameters(), fp_, edm::FileInPath::fullPath(), edm::EventSetup::get(), edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), estimatePileup::inputRange, edm::Service< T >::isAvailable(), j, minimumPosValue_, VarParsing::obj, printDebug_, printLog(), pset_, SiStripNoises::put(), matplotRender::reader, SiStripNoises::setData(), stripLengthMode_, and subDetAndLayer().

{
  // Read the gain from the given tag
  edm::ESHandle<SiStripApvGain> inputApvGain;
  iSetup.get<SiStripApvGainRcd>().get( inputApvGain );
  std::vector<uint32_t> inputDetIds;
  inputApvGain->getDetIds(inputDetIds);

  // Prepare the new object
  SiStripNoises* obj = new SiStripNoises();

  SiStripDetInfoFileReader reader(fp_.fullPath());



  stripLengthMode_ = pset_.getParameter<bool>("StripLengthMode");
  
  //parameters for random noise generation. not used if Strip length mode is chosen
  std::map<int, std::vector<double> > meanNoise;
  fillParameters(meanNoise, "MeanNoise");
  std::map<int, std::vector<double> > sigmaNoise;
  fillParameters(sigmaNoise, "SigmaNoise");
  minimumPosValue_ = pset_.getParameter<double>("MinPositiveNoise");

  //parameters for strip length proportional noise generation. not used if random mode is chosen
  std::map<int, std::vector<double> > noiseStripLengthLinearSlope;
  fillParameters(noiseStripLengthLinearSlope, "NoiseStripLengthSlope");
  std::map<int, std::vector<double> > noiseStripLengthLinearQuote;
  fillParameters(noiseStripLengthLinearQuote, "NoiseStripLengthQuote");
  electronsPerADC_ = pset_.getParameter<double>("electronPerAdc");        

  printDebug_ = pset_.getUntrackedParameter<uint32_t>("printDebug", 5);

  unsigned int count = 0;
  const std::map<uint32_t, SiStripDetInfoFileReader::DetInfo > DetInfos = reader.getAllData();
  for(std::map<uint32_t, SiStripDetInfoFileReader::DetInfo >::const_iterator it = DetInfos.begin(); it != DetInfos.end(); it++) {

    // Find if this DetId is in the input tag and if so how many are the Apvs for which it contains information
    SiStripApvGain::Range inputRange(inputApvGain->getRange(it->first));

    //Generate Noises for det detid
    SiStripNoises::InputVector theSiStripVector;
    float noise = 0.;
    uint32_t detId = it->first;
    std::pair<int, int> sl = subDetAndLayer(detId);
    unsigned short nApvs = it->second.nApvs;

    if(stripLengthMode_) {
      // Use strip length
      double linearSlope = noiseStripLengthLinearSlope[sl.first][sl.second];
      double linearQuote = noiseStripLengthLinearQuote[sl.first][sl.second];
      double stripLength = it->second.stripLength;
      for( unsigned short j=0; j<nApvs; ++j ) {

        double gain = inputApvGain->getApvGain(j, inputRange);

        for( unsigned short stripId = 0; stripId < 128; ++stripId ) {
          noise = ( ( linearSlope*stripLength + linearQuote) / electronsPerADC_ ) * gain;
          if( count<printDebug_ ) printLog(detId, stripId+128*j, noise);
          obj->setData(noise, theSiStripVector);
        }
      }
    }
    else {
      // Use random generator
      double meanN = meanNoise[sl.first][sl.second];
      double sigmaN = sigmaNoise[sl.first][sl.second];
      for( unsigned short j=0; j<nApvs; ++j ) {

        double gain = inputApvGain->getApvGain(j, inputRange);

        for( unsigned short stripId = 0; stripId < 128; ++stripId ) {
          noise = ( CLHEP::RandGauss::shoot(meanN, sigmaN) ) * gain;
          if( noise<=minimumPosValue_ ) noise = minimumPosValue_;
          if( count<printDebug_ ) printLog(detId, stripId+128*j, noise);
          obj->setData(noise, theSiStripVector);
        }
      }
    }
    ++count;
    
    if ( ! obj->put(it->first,theSiStripVector) ) {
      edm::LogError("SiStripNoisesFakeESSource::produce ")<<" detid already exists"<<std::endl;
    }
  }
  
  //End now write data in DB
  edm::Service<cond::service::PoolDBOutputService> mydbservice;
  
  if( mydbservice.isAvailable() ){
    if( mydbservice->isNewTagRequest("SiStripNoisesRcd") ){
      mydbservice->createNewIOV<SiStripNoises>(obj,mydbservice->beginOfTime(),mydbservice->endOfTime(),"SiStripNoisesRcd");
    }
    else {
      mydbservice->appendSinceTime<SiStripNoises>(obj,mydbservice->currentTime(),"SiStripNoisesRcd");
    }
  }
  else {
    edm::LogError("SiStripNoiseNormalizedWithApvGainBuilder")<<"Service is unavailable"<<std::endl;
  }
}
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 153 of file SiStripNoiseNormalizedWithApvGainBuilder.cc.

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

Referenced by analyze().

{
  int layersTIB = 4;
  int ringsTID = 3;
  int layersTOB = 6;
  int ringsTEC = 7;

  fillSubDetParameter( mapToFill, pset_.getParameter<std::vector<double> >(parameterName+"TIB"), int(StripSubdetector::TIB), layersTIB );
  fillSubDetParameter( mapToFill, pset_.getParameter<std::vector<double> >(parameterName+"TID"), int(StripSubdetector::TID), ringsTID );
  fillSubDetParameter( mapToFill, pset_.getParameter<std::vector<double> >(parameterName+"TOB"), int(StripSubdetector::TOB), layersTOB );
  fillSubDetParameter( mapToFill, pset_.getParameter<std::vector<double> >(parameterName+"TEC"), int(StripSubdetector::TEC), ringsTEC );
}
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 166 of file SiStripNoiseNormalizedWithApvGainBuilder.cc.

References Exception.

Referenced by fillParameters().

{
  if( v.size() == layers ) {
    mapToFill.insert(std::make_pair( subDet, v ));
  }
  else if( v.size() == 1 ) {
    std::vector<double> parV(layers, v[0]);
    mapToFill.insert(std::make_pair( subDet, parV ));
  }
  else {
    throw cms::Exception("Configuration") << "ERROR: number of parameters for subDet " << subDet << " are " << v.size() << ". They must be either 1 or " << layers << std::endl;
  }
}
void SiStripNoiseNormalizedWithApvGainBuilder::printLog ( const uint32_t  detId,
const unsigned short  strip,
const double &  noise 
) const [inline, private]

Definition at line 52 of file SiStripNoiseNormalizedWithApvGainBuilder.h.

Referenced by analyze().

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

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

Definition at line 127 of file SiStripNoiseNormalizedWithApvGainBuilder.cc.

References TIBDetId::layer(), TOBDetId::layer(), TIDDetId::ring(), TECDetId::ring(), DetId::subdetId(), StripSubdetector::TEC, StripSubdetector::TIB, StripSubdetector::TID, and StripSubdetector::TOB.

Referenced by analyze().

{
  int layerId = 0;

  StripSubdetector subid(detId);
  int subId = subid.subdetId();

  if( subId == int(StripSubdetector::TIB)) {
    TIBDetId theTIBDetId(detId);
    layerId = theTIBDetId.layer() - 1;
  }
  else if(subId == int(StripSubdetector::TOB)) {
    TOBDetId theTOBDetId(detId);
    layerId = theTOBDetId.layer() - 1;
  }
  else if(subId == int(StripSubdetector::TID)) {
    TIDDetId theTIDDetId(detId);
    layerId = theTIDDetId.ring() - 1;
  }
  if(subId == int(StripSubdetector::TEC)) {
    TECDetId theTECDetId = TECDetId(detId); 
    layerId = theTECDetId.ring() - 1;
  }
  return std::make_pair(subId, layerId);
}

Member Data Documentation

Definition at line 61 of file SiStripNoiseNormalizedWithApvGainBuilder.h.

Referenced by analyze().

Definition at line 57 of file SiStripNoiseNormalizedWithApvGainBuilder.h.

Referenced by analyze().

Definition at line 62 of file SiStripNoiseNormalizedWithApvGainBuilder.h.

Referenced by analyze().

Definition at line 58 of file SiStripNoiseNormalizedWithApvGainBuilder.h.

Definition at line 64 of file SiStripNoiseNormalizedWithApvGainBuilder.h.

Referenced by analyze().

Definition at line 59 of file SiStripNoiseNormalizedWithApvGainBuilder.h.

Referenced by analyze(), and fillParameters().

Definition at line 63 of file SiStripNoiseNormalizedWithApvGainBuilder.h.

Referenced by analyze().