CMS 3D CMS Logo

HPDNoiseLibraryReader Class Reference

#include <SimCalorimetry/HcalSimAlgos/interface/HPDNoiseLibraryReader.h>

List of all members.

Public Member Functions

std::vector< std::pair
< HcalDetId, const float * > > 
getBiasedNoisyHcalDetIds ()
std::vector< std::pair
< HcalDetId, const float * > > 
getNoisyHcalDetIds ()
 HPDNoiseLibraryReader (const edm::ParameterSet &)
 ~HPDNoiseLibraryReader ()

Static Public Member Functions

static void initializeServices ()

Public Attributes

std::string theHPDName
std::vector< std::string > theNames
std::vector< float > theNoiseRate
std::vector< inttheNoisyPhi
CLHEP::RandFlat * theRandFlat
CLHEP::RandGaussQ * theRandGaussQ
HPDNoiseReadertheReader
HcalTopology theTopology

Protected Member Functions

void setRandomEngine (CLHEP::HepRandomEngine &engine)
void setRandomEngine ()

Private Member Functions

bool applyNoise (int iphi)
void clearPhi ()
void fillRate ()
void getBiasedNoisyPhis ()
HPDNoiseDatagetNoiseData (int iphi)
void getNoisyPhis ()
std::string itos (int i)


Detailed Description

Definition at line 42 of file HPDNoiseLibraryReader.h.


Constructor & Destructor Documentation

HPDNoiseLibraryReader::HPDNoiseLibraryReader ( const edm::ParameterSet iConfig  ) 

Definition at line 22 of file HPDNoiseLibraryReader.cc.

References HPDNoiseReader::allNames(), fillRate(), edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), setRandomEngine(), theHPDName, theNames, and theReader.

00023 : theNoiseRate(0), theNoisyPhi(0), theRandFlat(0), theRandGaussQ(0)
00024 {
00025     setRandomEngine();
00026     
00027     ParameterSet pSet     = iConfig.getParameter<edm::ParameterSet>("HPDNoiseLibrary");
00028     FileInPath   filepath = pSet.getParameter<edm::FileInPath>("FileName");
00029     theHPDName            = pSet.getUntrackedParameter<string> ("HPDName", "HPD");
00030     string pName          = filepath.fullPath();
00031     if (pName.find(".") == 0) pName.erase(0,2);
00032     theReader = new HPDNoiseReader(pName);
00033     theNames = theReader->allNames();   // all 72x2 HPDs
00034 
00035     fillRate();
00036 }

HPDNoiseLibraryReader::~HPDNoiseLibraryReader (  ) 

Definition at line 38 of file HPDNoiseLibraryReader.cc.

References theRandFlat, and theRandGaussQ.

00039 {
00040     if (theRandFlat)delete theRandFlat;
00041     if (theRandGaussQ)delete theRandGaussQ;
00042 }


Member Function Documentation

bool HPDNoiseLibraryReader::applyNoise ( int  iphi  )  [private]

Definition at line 205 of file HPDNoiseLibraryReader.cc.

References find(), and theNoisyPhi.

Referenced by getNoiseData().

00206 {
00207 
00208     bool isAccepted = false;
00209     vector < int >::iterator phi_iter;
00210     phi_iter = find(theNoisyPhi.begin(), theNoisyPhi.end(), iphi);
00211     if (phi_iter != theNoisyPhi.end()){
00212         isAccepted = true;
00213     }
00214     return isAccepted;
00215 }

void HPDNoiseLibraryReader::clearPhi (  )  [private]

Definition at line 223 of file HPDNoiseLibraryReader.cc.

References theNoisyPhi.

Referenced by getBiasedNoisyPhis(), and getNoisyPhis().

00223                                      {
00224     theNoisyPhi.clear();
00225 }

void HPDNoiseLibraryReader::fillRate (  )  [private]

Definition at line 82 of file HPDNoiseLibraryReader.cc.

References TestMuL1L2Filter_cff::cerr, lat::endl(), HPDNoiseReader::getHandle(), i, HPDNoiseReader::rate(), theNames, theNoiseRate, theReader, and HPDNoiseReader::valid().

Referenced by HPDNoiseLibraryReader().

00083 {
00084     for (size_t i = 0; i < theNames.size(); ++i) {
00085         HPDNoiseReader::Handle hpdObj = theReader->getHandle(theNames[i]);
00086         if (theReader->valid(hpdObj)) {
00087             theNoiseRate.push_back(theReader->rate(hpdObj));
00088         } else {
00089             std::cerr << "HPD Handle Object is not valid!" << endl;
00090         }
00091     }
00092 }

vector< pair< HcalDetId, const float * > > HPDNoiseLibraryReader::getBiasedNoisyHcalDetIds (  ) 

Definition at line 184 of file HPDNoiseLibraryReader.cc.

References data, getBiasedNoisyPhis(), HPDNoiseData::getDataFrame(), HPDNoiseDataFrame::getFrame(), getNoiseData(), i, HPDNoiseDataFrame::id(), int, HLT_VtxMuL3::result, HPDNoiseData::size(), and theNoisyPhi.

00185 {
00186     
00187     vector< pair<HcalDetId, const float* > > result;
00188     // decide which phi are noisy by using noise rates 
00189     // and random numbers (to be called for each event)
00190     // at least one Phi is always noisy.
00191     getBiasedNoisyPhis();
00192     for (int i = 0; i < int(theNoisyPhi.size()); ++i) {
00193         int iphi = theNoisyPhi[i];
00194         HPDNoiseData* data;
00195         data = getNoiseData(iphi);
00196         for(unsigned int i=0; i<data->size();++i)
00197         {
00198             pair < HcalDetId, const float* >tmp_pair( data->getDataFrame(i).id(), data->getDataFrame(i).getFrame());
00199             result.push_back(tmp_pair);
00200         }
00201     }
00202     return result;
00203 }

void HPDNoiseLibraryReader::getBiasedNoisyPhis (  )  [private]

Definition at line 144 of file HPDNoiseLibraryReader.cc.

References clearPhi(), i, theNoiseRate, theNoisyPhi, and theRandFlat.

Referenced by getBiasedNoisyHcalDetIds().

00145 {
00146 
00147     clearPhi();
00148     double rndm[144];
00149 
00150     theRandFlat->shootArray(144, rndm);
00151     for (int i = 0; i < 144; ++i) {
00152         if (rndm[i] < theNoiseRate[i]) {
00153             theNoisyPhi.push_back(i + 1);
00154         }
00155     }
00156     // make sure one HPD is always noisy
00157     if (theNoisyPhi.size() == 0) {
00158         int iPhi = (theRandFlat->fireInt(144)) + 1; // integer from interval [0-144[ + 1
00159 
00160         theNoisyPhi.push_back(iPhi);
00161     }
00162 }

HPDNoiseData * HPDNoiseLibraryReader::getNoiseData ( int  iphi  )  [private]

Definition at line 94 of file HPDNoiseLibraryReader.cc.

References applyNoise(), TestMuL1L2Filter_cff::cerr, data, lat::endl(), HPDNoiseReader::getEntry(), HPDNoiseReader::getHandle(), itos(), name, theHPDName, theRandFlat, theReader, HPDNoiseReader::totalEntries(), and HPDNoiseReader::valid().

Referenced by getBiasedNoisyHcalDetIds(), and getNoisyHcalDetIds().

00095 {
00096     
00097     
00098     HPDNoiseData* data;
00099     // make sure that iphi from HcalDetId is found noisy at first.
00100     // In other words, be sure that iphi is in the collection of
00101     // noisy Phis
00102     if (!(applyNoise(iphi))) return data;
00103     
00104     int zside = 1;
00105     if(iphi>72){
00106        iphi = iphi-72;
00107        zside = -1;
00108     }
00109     std::string name;
00110     if (zside == 1) {
00111         name = "ZPlus" + theHPDName + itos(iphi);
00112     } else if (zside == -1){
00113         name = "ZMinus" + theHPDName + itos(iphi);
00114     }else {
00115         cerr << " ZSide Calculation Error." << endl;
00116     }
00117     HPDNoiseReader::Handle hpdObj = theReader->getHandle(name);
00118     if (theReader->valid(hpdObj)) {
00119         // randomly select one entry from library for this HPD
00120         unsigned long entry = theRandFlat->fireInt( theReader->totalEntries(hpdObj));
00121         theReader->getEntry(hpdObj, entry, &data);
00122     }else{
00123         std::cerr << " HPD Name in the library is not valid." << std::endl;
00124     }
00125     return data;
00126 }

vector< pair< HcalDetId, const float * > > HPDNoiseLibraryReader::getNoisyHcalDetIds (  ) 

Definition at line 164 of file HPDNoiseLibraryReader.cc.

References data, HPDNoiseData::getDataFrame(), HPDNoiseDataFrame::getFrame(), getNoiseData(), getNoisyPhis(), i, HPDNoiseDataFrame::id(), int, HLT_VtxMuL3::result, HPDNoiseData::size(), and theNoisyPhi.

00165 {
00166     
00167     vector< pair<HcalDetId, const float* > > result;
00168     // decide which phi are noisy by using noise rates 
00169     // and random numbers (to be called for each event)
00170     getNoisyPhis();
00171     for (int i = 0; i < int(theNoisyPhi.size()); ++i) {
00172         int iphi = theNoisyPhi[i];
00173         HPDNoiseData* data;
00174         data = getNoiseData(iphi);
00175         for(unsigned int i=0; i<data->size();++i)
00176         {
00177             pair < HcalDetId, const float* >tmp_pair( data->getDataFrame(i).id(), data->getDataFrame(i).getFrame());
00178             result.push_back(tmp_pair);
00179         }
00180     }
00181     return result;
00182 }

void HPDNoiseLibraryReader::getNoisyPhis (  )  [private]

Definition at line 129 of file HPDNoiseLibraryReader.cc.

References clearPhi(), i, theNoiseRate, theNoisyPhi, and theRandFlat.

Referenced by getNoisyHcalDetIds().

00130 {
00131 
00132     clearPhi();
00133     double rndm[144];
00134 
00135     theRandFlat->shootArray(144, rndm);
00136 
00137     for (int i = 0; i < 144; ++i) {
00138         if (rndm[i] < theNoiseRate[i]) {
00139             theNoisyPhi.push_back(i + 1);
00140         }
00141     }
00142 }

void HPDNoiseLibraryReader::initializeServices (  )  [static]

Definition at line 44 of file HPDNoiseLibraryReader.cc.

References edmplugin::standard::config(), edmplugin::PluginManager::configure(), edm::ServiceRegistry::createServicesFromConfig(), and edmplugin::PluginManager::isAvailable().

00045 {
00046     if (not edmplugin::PluginManager::isAvailable()) {
00047         edmplugin::PluginManager::configure(edmplugin::standard::config());
00048     }
00049 
00050     std::string config =
00051       "process CorrNoise = {"
00052           "service = RandomNumberGeneratorService" 
00053              "{" 
00054                   "untracked uint32 sourceSeed = 123456789" 
00055               "}" 
00056       "}";
00057 
00058     // create the services
00059     edm::ServiceToken tempToken = edm::ServiceRegistry::createServicesFromConfig(config);
00060 
00061     // make the services available
00062     edm::ServiceRegistry::Operate operate(tempToken);
00063 }

string HPDNoiseLibraryReader::itos ( int  i  )  [private]

Definition at line 216 of file HPDNoiseLibraryReader.cc.

References s.

Referenced by getNoiseData().

00217 {
00218     stringstream s;
00219     s << i;
00220     return s.str();
00221 }

void HPDNoiseLibraryReader::setRandomEngine ( CLHEP::HepRandomEngine &  engine  )  [protected]

Definition at line 74 of file HPDNoiseLibraryReader.cc.

References theRandFlat, and theRandGaussQ.

00075 {
00076     if(theRandGaussQ) delete theRandGaussQ;
00077     if(theRandFlat) delete theRandFlat;
00078     theRandGaussQ = new CLHEP::RandGaussQ(engine);
00079     theRandFlat = new CLHEP::RandFlat(engine);
00080 }

void HPDNoiseLibraryReader::setRandomEngine (  )  [protected]

Definition at line 64 of file HPDNoiseLibraryReader.cc.

References Exception, and edm::Service< T >::isAvailable().

Referenced by HPDNoiseLibraryReader().

00065 {
00066     edm::Service < edm::RandomNumberGenerator > rng;
00067     if (!rng.isAvailable()) {
00068         throw cms::Exception("Configuration") << "HcalHPDNoiseLibrary requires the RandomNumberGeneratorService\n"
00069           "which is not present in the configuration file.  You must add the service\n"
00070           "in the configuration file or remove the modules that require it.";
00071     }
00072     setRandomEngine(rng->getEngine());
00073 }


Member Data Documentation

std::string HPDNoiseLibraryReader::theHPDName

Definition at line 85 of file HPDNoiseLibraryReader.h.

Referenced by getNoiseData(), and HPDNoiseLibraryReader().

std::vector<std::string> HPDNoiseLibraryReader::theNames

Definition at line 84 of file HPDNoiseLibraryReader.h.

Referenced by fillRate(), and HPDNoiseLibraryReader().

std::vector<float> HPDNoiseLibraryReader::theNoiseRate

Definition at line 79 of file HPDNoiseLibraryReader.h.

Referenced by fillRate(), getBiasedNoisyPhis(), and getNoisyPhis().

std::vector<int> HPDNoiseLibraryReader::theNoisyPhi

Definition at line 80 of file HPDNoiseLibraryReader.h.

Referenced by applyNoise(), clearPhi(), getBiasedNoisyHcalDetIds(), getBiasedNoisyPhis(), getNoisyHcalDetIds(), and getNoisyPhis().

CLHEP::RandFlat* HPDNoiseLibraryReader::theRandFlat

Definition at line 81 of file HPDNoiseLibraryReader.h.

Referenced by getBiasedNoisyPhis(), getNoiseData(), getNoisyPhis(), setRandomEngine(), and ~HPDNoiseLibraryReader().

CLHEP::RandGaussQ* HPDNoiseLibraryReader::theRandGaussQ

Definition at line 82 of file HPDNoiseLibraryReader.h.

Referenced by setRandomEngine(), and ~HPDNoiseLibraryReader().

HPDNoiseReader* HPDNoiseLibraryReader::theReader

Definition at line 83 of file HPDNoiseLibraryReader.h.

Referenced by fillRate(), getNoiseData(), and HPDNoiseLibraryReader().

HcalTopology HPDNoiseLibraryReader::theTopology

Definition at line 76 of file HPDNoiseLibraryReader.h.


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:24:47 2009 for CMSSW by  doxygen 1.5.4