CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_13_patch3/src/SimCalorimetry/HcalSimAlgos/src/HPDNoiseLibraryReader.cc

Go to the documentation of this file.
00001 
00002 // --------------------------------------------------------
00003 // A class to read HPD noise from the library.
00004 // The deliverable of the class is the collection of
00005 // noisy HcalDetIds with associated noise in units of fC for
00006 // 10 time samples. During the library production a higher
00007 // theshold is used to find a noisy HPD. A lower threshold is
00008 // used to eliminate adding unnecessary quite channels to HPD 
00009 // noise event collection. Therefore user may not see whole 18 
00010 // channels for noisy HPD.
00011 //
00012 // Project: HPD noise library reader
00013 // Author: T.Yetkin University of Iowa, Feb. 7, 2008
00014 // $Id: HPDNoiseLibraryReader.cc,v 1.4 2011/06/17 12:34:22 abdullin Exp $
00015 // --------------------------------------------------------
00016 
00017 #include "SimCalorimetry/HcalSimAlgos/interface/HPDNoiseLibraryReader.h"
00018 
00019 using namespace edm;
00020 using namespace std;
00021 
00022 HPDNoiseLibraryReader::HPDNoiseLibraryReader(const edm::ParameterSet & iConfig)
00023 :theDischargeNoiseRate(0),
00024 theIonFeedbackFirstPeakRate(0),
00025 theIonFeedbackSecondPeakRate(0),
00026 theNoisyPhi(0),
00027 theRandFlat(0),
00028 theRandGaussQ(0) {
00029     setRandomEngine();
00030 
00031     ParameterSet pSet = iConfig.getParameter < edm::ParameterSet > ("HPDNoiseLibrary");
00032     FileInPath filepath = pSet.getParameter < edm::FileInPath > ("FileName");
00033 
00034     theHPDName = pSet.getUntrackedParameter < string > ("HPDName", "HPD");
00035     string pName = filepath.fullPath();
00036 
00037     if (pName.find(".") == 0)
00038         pName.erase(0, 2);
00039     theReader = new HPDNoiseReader(pName);
00040     theNames = theReader->allNames();   // all 72x2 HPDs
00041 
00042     fillRates();
00043 }
00044 
00045 HPDNoiseLibraryReader::~HPDNoiseLibraryReader() {
00046     if (theRandFlat)
00047         delete theRandFlat;
00048 
00049     if (theRandGaussQ)
00050         delete theRandGaussQ;
00051 }
00052 
00053 void HPDNoiseLibraryReader::initializeServices() {
00054     if (not edmplugin::PluginManager::isAvailable()) {
00055         edmplugin::PluginManager::configure(edmplugin::standard::config());
00056     }
00057 
00058     std::string config =
00059       "process CorrNoise = {"
00060       "service = RandomNumberGeneratorService" "{" "untracked uint32 sourceSeed = 123456789" "}" "}";
00061 
00062     // create the services
00063     edm::ServiceToken tempToken = edm::ServiceRegistry::createServicesFromConfig(config);
00064 
00065     // make the services available
00066     edm::ServiceRegistry::Operate operate(tempToken);
00067 }
00068 void HPDNoiseLibraryReader::setRandomEngine() {
00069     edm::Service < edm::RandomNumberGenerator > rng;
00070     if (!rng.isAvailable()) {
00071         throw cms::Exception("Configuration") << "HcalHPDNoiseLibrary requires the RandomNumberGeneratorService\n"
00072           "which is not present in the configuration file.  You must add the service\n"
00073           "in the configuration file or remove the modules that require it.";
00074     }
00075     setRandomEngine(rng->getEngine());
00076 }
00077 void HPDNoiseLibraryReader::setRandomEngine(CLHEP::HepRandomEngine & engine) {
00078     if (theRandGaussQ)
00079         delete theRandGaussQ;
00080 
00081     if (theRandFlat)
00082         delete theRandFlat;
00083 
00084     theRandGaussQ = new CLHEP::RandGaussQ(engine);
00085     theRandFlat = new CLHEP::RandFlat(engine);
00086 }
00087 
00088 void HPDNoiseLibraryReader::fillRates() {
00089     for (size_t i = 0; i < theNames.size(); ++i) {
00090         HPDNoiseReader::Handle hpdObj = theReader->getHandle(theNames[i]);
00091         if (theReader->valid(hpdObj)) {
00092             theDischargeNoiseRate.push_back(theReader->dischargeRate(hpdObj));
00093             theIonFeedbackFirstPeakRate.push_back(theReader->ionFeedbackFirstPeakRate(hpdObj));
00094             theIonFeedbackSecondPeakRate.push_back(theReader->ionFeedbackSecondPeakRate(hpdObj));
00095         } else {
00096             std::cerr << "HPD Handle Object is not valid!" << endl;
00097         }
00098     }
00099 }
00100 
00101 HPDNoiseData *HPDNoiseLibraryReader::getNoiseData(int iphi) {
00102 
00103 
00104     HPDNoiseData *data = 0; //data->size() is checked wherever actually used  
00105 
00106     // make sure that iphi from HcalDetId is found noisy at first.
00107     // In other words, be sure that iphi is in the collection of
00108     // noisy Phis
00109     if (!(IsNoiseApplicable(iphi)))         
00110         return data;
00111 
00112     int zside = 1;
00113 
00114     if (iphi > 72) {
00115         iphi = iphi - 72;
00116         zside = -1;
00117     }
00118     std::string name;
00119     if (zside == 1) {
00120         name = "ZPlus" + theHPDName + itos(iphi);
00121     } else if (zside == -1) {
00122         name = "ZMinus" + theHPDName + itos(iphi);
00123     } else {
00124         cerr << " ZSide Calculation Error." << endl;
00125     }
00126     HPDNoiseReader::Handle hpdObj = theReader->getHandle(name);
00127     if (theReader->valid(hpdObj)) {
00128         // randomly select one entry from library for this HPD
00129         unsigned long entry = theRandFlat->fireInt(theReader->totalEntries(hpdObj));
00130 
00131         theReader->getEntry(hpdObj, entry, &data);
00132     } else {
00133         std::cerr << " HPD Name in the library is not valid." << std::endl;
00134     }
00135     return data;
00136 }
00137 
00138 
00139 void HPDNoiseLibraryReader::getNoisyPhis() {
00140 
00141     clearPhi();
00142     double rndm[144];
00143 
00144     theRandFlat->shootArray(144, rndm);
00145 
00146     for (int i = 0; i < 144; ++i) {
00147         if (rndm[i] < theDischargeNoiseRate[i]) {
00148             theNoisyPhi.push_back(i + 1);
00149         }
00150     }
00151 }
00152 
00153 void HPDNoiseLibraryReader::getBiasedNoisyPhis() {
00154 
00155     clearPhi();
00156     double rndm[144];
00157 
00158     theRandFlat->shootArray(144, rndm);
00159     for (int i = 0; i < 144; ++i) {
00160         if (rndm[i] < theDischargeNoiseRate[i]) {
00161             theNoisyPhi.push_back(i + 1);
00162         }
00163     }
00164     // make sure one HPD is always noisy
00165     if (theNoisyPhi.size() == 0) {
00166         int iPhi = (theRandFlat->fireInt(144)) + 1; // integer from interval [0-144[ + 1
00167 
00168         theNoisyPhi.push_back(iPhi);
00169     }
00170 }
00171 
00172 vector <pair<HcalDetId, const float *> >HPDNoiseLibraryReader::getNoisyHcalDetIds() {
00173 
00174     vector <pair< HcalDetId, const float *> >result;
00175 
00176     // decide which phi are noisy by using noise rates 
00177     // and random numbers (to be called for each event)
00178     getNoisyPhis();
00179     for (int i = 0; i < int (theNoisyPhi.size()); ++i) {
00180         int iphi = theNoisyPhi[i];
00181         HPDNoiseData *data;
00182 
00183         data = getNoiseData(iphi);
00184         for (unsigned int i = 0; i < data->size(); ++i) {
00185             pair < HcalDetId, const float *>tmp_pair(data->getDataFrame(i).id(), data->getDataFrame(i).getFrame());
00186 
00187             result.push_back(tmp_pair);
00188         }
00189     }
00190     return result;
00191 }
00192 
00193 vector <pair<HcalDetId, const float *> >HPDNoiseLibraryReader::getNoisyHcalDetIds(int timeSliceId) 
00194 {
00195     vector <pair< HcalDetId, const float *> >result;
00196     // decide which phi are noisy by using noise rates 
00197     // and random numbers (to be called for each event)
00198     getNoisyPhis();
00199     for (int i = 0; i < int (theNoisyPhi.size()); ++i) {
00200         int iphi = theNoisyPhi[i];
00201         HPDNoiseData *data;
00202 
00203         data = getNoiseData(iphi);
00204         for (unsigned int i = 0; i < data->size(); ++i) {
00205             float* data_ = const_cast<float*>(data->getDataFrame(i).getFrame());
00206             shuffleData(timeSliceId, data_);
00207             const float* _data_ =const_cast<const float*>(data_);
00208             pair < HcalDetId, const float *>tmp_pair(data->getDataFrame(i).id(), _data_);
00209             result.push_back(tmp_pair);
00210         }
00211     }
00212     return result;
00213 
00214 }
00215 vector < pair < HcalDetId, const float *> >HPDNoiseLibraryReader::getBiasedNoisyHcalDetIds(int timeSliceId) {
00216 
00217     vector < pair < HcalDetId, const float *> >result;
00218 
00219     // decide which phi are noisy by using noise rates 
00220     // and random numbers (to be called for each event)
00221     // at least one Phi is always noisy.
00222     getBiasedNoisyPhis();
00223     for (int i = 0; i < int (theNoisyPhi.size()); ++i) {
00224         int iphi = theNoisyPhi[i];
00225         HPDNoiseData *data;
00226 
00227         data = getNoiseData(iphi);
00228         for (unsigned int i = 0; i < data->size(); ++i) {
00229             float* data_ = const_cast<float*>(data->getDataFrame(i).getFrame());
00230             shuffleData(timeSliceId, data_);
00231             const float* _data_ =const_cast<const float*>(data_);
00232             pair < HcalDetId, const float *>tmp_pair(data->getDataFrame(i).id(), _data_);
00233             result.push_back(tmp_pair);
00234         }
00235     }
00236     return result;
00237 }
00238 
00239 vector < pair < HcalDetId, const float *> >HPDNoiseLibraryReader::getBiasedNoisyHcalDetIds() {
00240 
00241     vector < pair < HcalDetId, const float *> >result;
00242 
00243     // decide which phi are noisy by using noise rates 
00244     // and random numbers (to be called for each event)
00245     // at least one Phi is always noisy.
00246     getBiasedNoisyPhis();
00247     for (int i = 0; i < int (theNoisyPhi.size()); ++i) {
00248         int iphi = theNoisyPhi[i];
00249         HPDNoiseData *data;
00250 
00251         data = getNoiseData(iphi);
00252         for (unsigned int i = 0; i < data->size(); ++i) {
00253             pair < HcalDetId, const float *>tmp_pair(data->getDataFrame(i).id(), data->getDataFrame(i).getFrame());
00254 
00255             result.push_back(tmp_pair);
00256         }
00257     }
00258     return result;
00259 }
00260 
00261 double HPDNoiseLibraryReader::getIonFeedbackNoise(HcalDetId id, double energy, double bias) {
00262 
00263     // constants for simulation/parameterization
00264     double pe2Charge = 0.333333;    // fC/p.e.
00265     double GeVperfC = 0.177;    // in unit GeV/fC and this will be updated when it start reading from DB.
00266     double PedSigma = 0.8;
00267     double noise = 0.;          // fC
00268 
00269     int iphi = (id.ieta() > 0) ? (id.iphi()) : (id.iphi() + 72);
00270     double rateInTail = theIonFeedbackFirstPeakRate[iphi - 1];
00271     double rateInSecondTail = theIonFeedbackSecondPeakRate[iphi - 1];
00272 
00273     if (bias != 0.) {
00274         rateInTail = rateInTail * bias;
00275         rateInSecondTail = rateInSecondTail * bias;
00276     } else {
00277         edm::LogError("HPDNoise") << "HPDNoise: ion feedback error (biased or unbiased selection)." << bias << " failed";
00278         throw cms::Exception("Unknown", "biase selection ")
00279         << "Usage of " << bias << " fails\n";
00280     }
00281     double Charge = energy / GeVperfC;
00282 
00283     // three gauss fit is applied to data to get ion feedback distribution
00284     // the script is at neutralino: /home/tyetkin/work/hpd_noise/PlotFromPelin.C
00285     // a new fit woth double-sigmoids is under way.
00286     // parameters (in fC)
00287     // first gaussian
00288     // double p0 = 9.53192e+05;
00289     // double p1 = -3.13653e-01;
00290     // double p2 = 2.78350e+00;
00291 
00292     // second gaussian
00293     // double p3 = 2.41611e+03;
00294     double p4 = 2.06117e+01;
00295     double p5 = 1.09239e+01;
00296 
00297     // third gaussian
00298     // double p6 = 3.42793e+01;
00299     double p7 = 5.45548e+01;
00300     double p8 = 1.59696e+01;
00301 
00302     if (Charge > 3 * PedSigma) {    // 3 sigma away from pedestal mean
00303         int npe = int (Charge / pe2Charge);
00304         double a = 0.;
00305         double b = 0.;
00306 
00307         for (int j = 0; j < npe; ++j) {
00308             double probability = theRandFlat->shoot();
00309 
00310             if (probability < rateInTail) { // total tail
00311                 if (probability < rateInSecondTail) {   // second tail
00312                     Rannor(a, b);
00313                     noise += b * p8 + p7;
00314                 } else {
00315                     Rannor(a, b);
00316                     noise += b * p5 + p4;
00317                 }
00318             }
00319         }
00320         // add pedestal 
00321         if (noise > 0.)
00322             noise += theRandGaussQ->fire(0, 2 * PedSigma);
00323     }
00324     return (noise * GeVperfC);  // returns noise in GeV.
00325 
00326 }
00327 
00328 bool HPDNoiseLibraryReader::IsNoiseApplicable(int iphi) {
00329 
00330     bool isAccepted = false;
00331     vector < int >::iterator phi_iter;
00332 
00333     phi_iter = find(theNoisyPhi.begin(), theNoisyPhi.end(), iphi);
00334     if (phi_iter != theNoisyPhi.end()) {
00335         isAccepted = true;
00336     }
00337     return isAccepted;
00338 }
00339 void HPDNoiseLibraryReader::shuffleData(int timeSliceId, float* &data)
00340 {
00341    if(timeSliceId == -1 || (timeSliceId>=10)) return;
00342    //make a local copy of input data
00343    float Data[10] = {0.,0.,0.,0.,0.,0.,0.,0.,0.,0.};
00344    for(int i=0;i<10;++i){
00345        Data[i] = data[i];
00346    }
00347    int ts_max = -1;
00348    float max = -999.;
00349    for(int i=0;i<10;++i){
00350        if(Data[i]>max){
00351            max = data[i];
00352            ts_max = i;
00353        }
00354    }
00355    if((ts_max == -1)){//couldn't find ts_max, return the same value.
00356        return;
00357    }else{
00358        // always shift the noise to the right by putting zeroes to the previous slices.
00359        // the noise is pedestal subtracted. 0 value is acceptable.
00360        int k = -1;
00361        for(int i=0;i<10;++i){
00362            data[i] = 0.;
00363            int newIdx = timeSliceId+k;
00364            float dd = 0.;
00365            if(newIdx < 10){
00366                data[newIdx] = Data[ts_max+k];
00367                dd = Data[ts_max+k];
00368                i = newIdx;
00369            }
00370            data[i] = dd;
00371            ++k;
00372        }
00373                                                                                                            
00374    }
00375 }
00376 
00377 //I couldn't find Rannor in CLHEP/Random. For now, use it from ROOT (copy/paste) by little modification.
00378 void HPDNoiseLibraryReader::Rannor(double &a, double &b) {
00379     double r,
00380       x,
00381       y,
00382       z;
00383 
00384     y = theRandFlat->shoot();
00385     z = theRandFlat->shoot();
00386     x = z * 6.28318530717958623;
00387     r = TMath::Sqrt(-2 * TMath::Log(y));
00388     a = r * TMath::Sin(x);
00389     b = r * TMath::Cos(x);
00390 }
00391 string HPDNoiseLibraryReader::itos(int i) {
00392     stringstream s;
00393 
00394     s << i;
00395     return s.str();
00396 }
00397 
00398 void HPDNoiseLibraryReader::clearPhi() {
00399     theNoisyPhi.clear();
00400 }