CMS 3D CMS Logo

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

#include <HPDNoiseLibraryReader.h>

Public Member Functions

std::vector< std::pair
< HcalDetId, const float * > > 
getBiasedNoisyHcalDetIds ()
 
std::vector< std::pair
< HcalDetId, const float * > > 
getBiasedNoisyHcalDetIds (int timeSliceId)
 
double getIonFeedbackNoise (HcalDetId id, double energy, double bias)
 
std::vector< std::pair
< HcalDetId, const float * > > 
getNoisyHcalDetIds ()
 
std::vector< std::pair
< HcalDetId, const float * > > 
getNoisyHcalDetIds (int timeSliceId)
 
 HPDNoiseLibraryReader (const edm::ParameterSet &)
 
 ~HPDNoiseLibraryReader ()
 

Static Public Member Functions

static void initializeServices ()
 

Public Attributes

std::vector< float > theDischargeNoiseRate
 
std::string theHPDName
 
std::vector< float > theIonFeedbackFirstPeakRate
 
std::vector< float > theIonFeedbackSecondPeakRate
 
std::vector< std::string > theNames
 
std::vector< int > theNoisyPhi
 
CLHEP::RandFlat * theRandFlat
 
CLHEP::RandGaussQ * theRandGaussQ
 
HPDNoiseReadertheReader
 

Protected Member Functions

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

Private Member Functions

void clearPhi ()
 
void fillRates ()
 
void getBiasedNoisyPhis ()
 
HPDNoiseDatagetNoiseData (int iphi)
 
void getNoisyPhis ()
 
bool IsNoiseApplicable (int iphi)
 
std::string itos (int i)
 
void Rannor (double &a, double &b)
 
void shuffleData (int timeSliceId, float *&data)
 

Detailed Description

Definition at line 43 of file HPDNoiseLibraryReader.h.

Constructor & Destructor Documentation

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

Definition at line 21 of file HPDNoiseLibraryReader.cc.

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

25 theNoisyPhi(0),
26 theRandFlat(0),
27 theRandGaussQ(0) {
29 
30  ParameterSet pSet = iConfig.getParameter < edm::ParameterSet > ("HPDNoiseLibrary");
31  FileInPath filepath = pSet.getParameter < edm::FileInPath > ("FileName");
32 
33  theHPDName = pSet.getUntrackedParameter < string > ("HPDName", "HPD");
34  string pName = filepath.fullPath();
35 
36  if (pName.find(".") == 0)
37  pName.erase(0, 2);
38  theReader = new HPDNoiseReader(pName);
39  theNames = theReader->allNames(); // all 72x2 HPDs
40 
41  fillRates();
42 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
std::vector< float > theIonFeedbackFirstPeakRate
CLHEP::RandGaussQ * theRandGaussQ
std::vector< int > theNoisyPhi
std::vector< float > theIonFeedbackSecondPeakRate
std::vector< std::string > theNames
std::vector< std::string > allNames() const
all HPD instances in the library
std::vector< float > theDischargeNoiseRate
CLHEP::RandFlat * theRandFlat
HPDNoiseLibraryReader::~HPDNoiseLibraryReader ( )

Definition at line 44 of file HPDNoiseLibraryReader.cc.

References theRandFlat, and theRandGaussQ.

44  {
45  if (theRandFlat)
46  delete theRandFlat;
47 
48  if (theRandGaussQ)
49  delete theRandGaussQ;
50 }
CLHEP::RandGaussQ * theRandGaussQ
CLHEP::RandFlat * theRandFlat

Member Function Documentation

void HPDNoiseLibraryReader::clearPhi ( )
private

Definition at line 392 of file HPDNoiseLibraryReader.cc.

References theNoisyPhi.

Referenced by getBiasedNoisyPhis(), and getNoisyPhis().

392  {
393  theNoisyPhi.clear();
394 }
std::vector< int > theNoisyPhi
void HPDNoiseLibraryReader::fillRates ( )
private

Definition at line 87 of file HPDNoiseLibraryReader.cc.

References dtNoiseDBValidation_cfg::cerr, HPDNoiseReader::dischargeRate(), HPDNoiseReader::getHandle(), i, HPDNoiseReader::ionFeedbackFirstPeakRate(), HPDNoiseReader::ionFeedbackSecondPeakRate(), theDischargeNoiseRate, theIonFeedbackFirstPeakRate, theIonFeedbackSecondPeakRate, theNames, theReader, and HPDNoiseReader::valid().

Referenced by HPDNoiseLibraryReader().

87  {
88  for (size_t i = 0; i < theNames.size(); ++i) {
90  if (theReader->valid(hpdObj)) {
91  theDischargeNoiseRate.push_back(theReader->dischargeRate(hpdObj));
94  } else {
95  std::cerr << "HPD Handle Object is not valid!" << endl;
96  }
97  }
98 }
int i
Definition: DBlmapReader.cc:9
std::vector< float > theIonFeedbackFirstPeakRate
float dischargeRate(Handle fHandle) const
discharge rate for the instance
Handle getHandle(const std::string &fName)
get handle to access data for one HPD instance
float ionFeedbackSecondPeakRate(Handle fHandle) const
std::vector< float > theIonFeedbackSecondPeakRate
bool valid(Handle fHandle) const
test if handle is valid
std::vector< std::string > theNames
std::vector< float > theDischargeNoiseRate
float ionFeedbackFirstPeakRate(Handle fHandle) const
ionfeedback rate for the instance
vector< pair< HcalDetId, const float * > > HPDNoiseLibraryReader::getBiasedNoisyHcalDetIds ( )

Definition at line 235 of file HPDNoiseLibraryReader.cc.

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

235  {
236 
237  vector < pair < HcalDetId, const float *> >result;
238 
239  // decide which phi are noisy by using noise rates
240  // and random numbers (to be called for each event)
241  // at least one Phi is always noisy.
243  for (int i = 0; i < int (theNoisyPhi.size()); ++i) {
244  int iphi = theNoisyPhi[i];
246 
247  data = getNoiseData(iphi);
248  for (unsigned int i = 0; i < data->size(); ++i) {
249  result.emplace_back(data->getDataFrame(i).id(), data->getDataFrame(i).getFrame());
250  }
251  }
252  return result;
253 }
int i
Definition: DBlmapReader.cc:9
const HPDNoiseDataFrame & getDataFrame(size_t i) const
retrive frame for the given index
Definition: HPDNoiseData.cc:21
const float * getFrame() const
array of 10 charges corresponding to one channel
std::vector< int > theNoisyPhi
tuple result
Definition: query.py:137
HcalDetId id() const
detId for the frame
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
HPDNoiseData * getNoiseData(int iphi)
unsigned size() const
number of noise channels in the event
Definition: HPDNoiseData.h:24
vector< pair< HcalDetId, const float * > > HPDNoiseLibraryReader::getBiasedNoisyHcalDetIds ( int  timeSliceId)

Definition at line 212 of file HPDNoiseLibraryReader.cc.

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

212  {
213 
214  vector < pair < HcalDetId, const float *> >result;
215 
216  // decide which phi are noisy by using noise rates
217  // and random numbers (to be called for each event)
218  // at least one Phi is always noisy.
220  for (int i = 0; i < int (theNoisyPhi.size()); ++i) {
221  int iphi = theNoisyPhi[i];
223 
224  data = getNoiseData(iphi);
225  for (unsigned int i = 0; i < data->size(); ++i) {
226  float* data_ = const_cast<float*>(data->getDataFrame(i).getFrame());
227  shuffleData(timeSliceId, data_);
228  const float* _data_ =const_cast<const float*>(data_);
229  result.emplace_back(data->getDataFrame(i).id(), _data_);
230  }
231  }
232  return result;
233 }
int i
Definition: DBlmapReader.cc:9
const HPDNoiseDataFrame & getDataFrame(size_t i) const
retrive frame for the given index
Definition: HPDNoiseData.cc:21
const float * getFrame() const
array of 10 charges corresponding to one channel
std::vector< int > theNoisyPhi
tuple result
Definition: query.py:137
HcalDetId id() const
detId for the frame
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
HPDNoiseData * getNoiseData(int iphi)
unsigned size() const
number of noise channels in the event
Definition: HPDNoiseData.h:24
void shuffleData(int timeSliceId, float *&data)
void HPDNoiseLibraryReader::getBiasedNoisyPhis ( )
private

Definition at line 152 of file HPDNoiseLibraryReader.cc.

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

Referenced by getBiasedNoisyHcalDetIds().

152  {
153 
154  clearPhi();
155  double rndm[144];
156 
157  theRandFlat->shootArray(144, rndm);
158  for (int i = 0; i < 144; ++i) {
159  if (rndm[i] < theDischargeNoiseRate[i]) {
160  theNoisyPhi.push_back(i + 1);
161  }
162  }
163  // make sure one HPD is always noisy
164  if (theNoisyPhi.size() == 0) {
165  int iPhi = (theRandFlat->fireInt(144)) + 1; // integer from interval [0-144[ + 1
166 
167  theNoisyPhi.push_back(iPhi);
168  }
169 }
int i
Definition: DBlmapReader.cc:9
std::vector< int > theNoisyPhi
std::vector< float > theDischargeNoiseRate
CLHEP::RandFlat * theRandFlat
double HPDNoiseLibraryReader::getIonFeedbackNoise ( HcalDetId  id,
double  energy,
double  bias 
)

HPD Ion feedback simulation based on LED data. A simple simulation which uses gaussian fit to data. biased = false ==> HPD noise from Ion Feedback only, unbiased biased = true ==> HPD noise from Ion Feedback only, biased (rate is X times larger than nominal rate)

Definition at line 255 of file HPDNoiseLibraryReader.cc.

References a, b, edm::hlt::Exception, j, p4, pe2Charge, Rannor(), theIonFeedbackFirstPeakRate, theIonFeedbackSecondPeakRate, theRandFlat, and theRandGaussQ.

255  {
256 
257  // constants for simulation/parameterization
258  double pe2Charge = 0.333333; // fC/p.e.
259  double GeVperfC = 0.177; // in unit GeV/fC and this will be updated when it start reading from DB.
260  double PedSigma = 0.8;
261  double noise = 0.; // fC
262 
263  int iphi = (id.ieta() > 0) ? (id.iphi()) : (id.iphi() + 72);
264  double rateInTail = theIonFeedbackFirstPeakRate[iphi - 1];
265  double rateInSecondTail = theIonFeedbackSecondPeakRate[iphi - 1];
266 
267  if (bias != 0.) {
268  rateInTail = rateInTail * bias;
269  rateInSecondTail = rateInSecondTail * bias;
270  } else {
271  edm::LogError("HPDNoise") << "HPDNoise: ion feedback error (biased or unbiased selection)." << bias << " failed";
272  throw cms::Exception("Unknown", "biase selection ")
273  << "Usage of " << bias << " fails\n";
274  }
275  double Charge = energy / GeVperfC;
276 
277  // three gauss fit is applied to data to get ion feedback distribution
278  // the script is at neutralino: /home/tyetkin/work/hpd_noise/PlotFromPelin.C
279  // a new fit woth double-sigmoids is under way.
280  // parameters (in fC)
281  // first gaussian
282  // double p0 = 9.53192e+05;
283  // double p1 = -3.13653e-01;
284  // double p2 = 2.78350e+00;
285 
286  // second gaussian
287  // double p3 = 2.41611e+03;
288  double p4 = 2.06117e+01;
289  double p5 = 1.09239e+01;
290 
291  // third gaussian
292  // double p6 = 3.42793e+01;
293  double p7 = 5.45548e+01;
294  double p8 = 1.59696e+01;
295 
296  if (Charge > 3 * PedSigma) { // 3 sigma away from pedestal mean
297  int npe = int (Charge / pe2Charge);
298  double a = 0.;
299  double b = 0.;
300 
301  for (int j = 0; j < npe; ++j) {
302  double probability = theRandFlat->shoot();
303 
304  if (probability < rateInTail) { // total tail
305  if (probability < rateInSecondTail) { // second tail
306  Rannor(a, b);
307  noise += b * p8 + p7;
308  } else {
309  Rannor(a, b);
310  noise += b * p5 + p4;
311  }
312  }
313  }
314  // add pedestal
315  if (noise > 0.)
316  noise += theRandGaussQ->fire(0, 2 * PedSigma);
317  }
318  return (noise * GeVperfC); // returns noise in GeV.
319 
320 }
std::vector< float > theIonFeedbackFirstPeakRate
CLHEP::RandGaussQ * theRandGaussQ
double pe2Charge
std::vector< float > theIonFeedbackSecondPeakRate
double p4[4]
Definition: TauolaWrapper.h:92
int j
Definition: DBlmapReader.cc:9
void Rannor(double &a, double &b)
double b
Definition: hdecay.h:120
CLHEP::RandFlat * theRandFlat
double a
Definition: hdecay.h:121
HPDNoiseData * HPDNoiseLibraryReader::getNoiseData ( int  iphi)
private

Definition at line 100 of file HPDNoiseLibraryReader.cc.

References dtNoiseDBValidation_cfg::cerr, data, HPDNoiseReader::getEntry(), HPDNoiseReader::getHandle(), IsNoiseApplicable(), itos(), mergeVDriftHistosByStation::name, AlCaHLTBitMon_QueryRunRegistry::string, theHPDName, theRandFlat, theReader, HPDNoiseReader::totalEntries(), and HPDNoiseReader::valid().

Referenced by getBiasedNoisyHcalDetIds(), and getNoisyHcalDetIds().

100  {
101 
102 
103  HPDNoiseData *data = 0; //data->size() is checked wherever actually used
104 
105  // make sure that iphi from HcalDetId is found noisy at first.
106  // In other words, be sure that iphi is in the collection of
107  // noisy Phis
108  if (!(IsNoiseApplicable(iphi)))
109  return data;
110 
111  int zside = 1;
112 
113  if (iphi > 72) {
114  iphi = iphi - 72;
115  zside = -1;
116  }
118  if (zside == 1) {
119  name = "ZPlus" + theHPDName + itos(iphi);
120  } else if (zside == -1) {
121  name = "ZMinus" + theHPDName + itos(iphi);
122  } else {
123  cerr << " ZSide Calculation Error." << endl;
124  }
126  if (theReader->valid(hpdObj)) {
127  // randomly select one entry from library for this HPD
128  unsigned long entry = theRandFlat->fireInt(theReader->totalEntries(hpdObj));
129 
130  theReader->getEntry(hpdObj, entry, &data);
131  } else {
132  std::cerr << " HPD Name in the library is not valid." << std::endl;
133  }
134  return data;
135 }
Handle getHandle(const std::string &fName)
get handle to access data for one HPD instance
unsigned long totalEntries(Handle fHandle) const
total number of noise events for the HPD instance in the library
bool valid(Handle fHandle) const
test if handle is valid
CLHEP::RandFlat * theRandFlat
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
void getEntry(Handle fHandle, HPDNoiseData **fData)
retrive one entry from the sequentially
vector< pair< HcalDetId, const float * > > HPDNoiseLibraryReader::getNoisyHcalDetIds ( )

Definition at line 171 of file HPDNoiseLibraryReader.cc.

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

Referenced by HPDNoiseGenerator::fillNoiseSignals().

171  {
172 
173  vector <pair< HcalDetId, const float *> >result;
174 
175  // decide which phi are noisy by using noise rates
176  // and random numbers (to be called for each event)
177  getNoisyPhis();
178  for (int i = 0; i < int (theNoisyPhi.size()); ++i) {
179  int iphi = theNoisyPhi[i];
181 
182  data = getNoiseData(iphi);
183  for (unsigned int i = 0; i < data->size(); ++i) {
184 
185  result.emplace_back(data->getDataFrame(i).id(), data->getDataFrame(i).getFrame());
186  }
187  }
188  return result;
189 }
int i
Definition: DBlmapReader.cc:9
const HPDNoiseDataFrame & getDataFrame(size_t i) const
retrive frame for the given index
Definition: HPDNoiseData.cc:21
const float * getFrame() const
array of 10 charges corresponding to one channel
std::vector< int > theNoisyPhi
tuple result
Definition: query.py:137
HcalDetId id() const
detId for the frame
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
HPDNoiseData * getNoiseData(int iphi)
unsigned size() const
number of noise channels in the event
Definition: HPDNoiseData.h:24
vector< pair< HcalDetId, const float * > > HPDNoiseLibraryReader::getNoisyHcalDetIds ( int  timeSliceId)

Definition at line 191 of file HPDNoiseLibraryReader.cc.

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

192 {
193  vector <pair< HcalDetId, const float *> >result;
194  // decide which phi are noisy by using noise rates
195  // and random numbers (to be called for each event)
196  getNoisyPhis();
197  for (int i = 0; i < int (theNoisyPhi.size()); ++i) {
198  int iphi = theNoisyPhi[i];
200 
201  data = getNoiseData(iphi);
202  for (unsigned int i = 0; i < data->size(); ++i) {
203  float* data_ = const_cast<float*>(data->getDataFrame(i).getFrame());
204  shuffleData(timeSliceId, data_);
205  const float* _data_ =const_cast<const float*>(data_);
206  result.emplace_back(data->getDataFrame(i).id(), _data_);
207  }
208  }
209  return result;
210 
211 }
int i
Definition: DBlmapReader.cc:9
const HPDNoiseDataFrame & getDataFrame(size_t i) const
retrive frame for the given index
Definition: HPDNoiseData.cc:21
const float * getFrame() const
array of 10 charges corresponding to one channel
std::vector< int > theNoisyPhi
tuple result
Definition: query.py:137
HcalDetId id() const
detId for the frame
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
HPDNoiseData * getNoiseData(int iphi)
unsigned size() const
number of noise channels in the event
Definition: HPDNoiseData.h:24
void shuffleData(int timeSliceId, float *&data)
void HPDNoiseLibraryReader::getNoisyPhis ( )
private

Definition at line 138 of file HPDNoiseLibraryReader.cc.

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

Referenced by getNoisyHcalDetIds().

138  {
139 
140  clearPhi();
141  double rndm[144];
142 
143  theRandFlat->shootArray(144, rndm);
144 
145  for (int i = 0; i < 144; ++i) {
146  if (rndm[i] < theDischargeNoiseRate[i]) {
147  theNoisyPhi.push_back(i + 1);
148  }
149  }
150 }
int i
Definition: DBlmapReader.cc:9
std::vector< int > theNoisyPhi
std::vector< float > theDischargeNoiseRate
CLHEP::RandFlat * theRandFlat
void HPDNoiseLibraryReader::initializeServices ( )
static

Definition at line 52 of file HPDNoiseLibraryReader.cc.

References HDQMDatabaseProducer::config, edmplugin::standard::config(), edmplugin::PluginManager::configure(), edm::ServiceRegistry::createServicesFromConfig(), edmplugin::PluginManager::isAvailable(), cmsPerfStripChart::operate(), and AlCaHLTBitMon_QueryRunRegistry::string.

52  {
55  }
56 
58  "process CorrNoise = {"
59  "service = RandomNumberGeneratorService" "{" "untracked uint32 sourceSeed = 123456789" "}" "}";
60 
61  // create the services
63 
64  // make the services available
66 }
static PluginManager & configure(const Config &)
static ServiceToken createServicesFromConfig(std::string const &config)
PluginManager::Config config()
Definition: standard.cc:21
bool HPDNoiseLibraryReader::IsNoiseApplicable ( int  iphi)
private

Definition at line 322 of file HPDNoiseLibraryReader.cc.

References spr::find(), and theNoisyPhi.

Referenced by getNoiseData().

322  {
323 
324  bool isAccepted = false;
325  vector < int >::iterator phi_iter;
326 
327  phi_iter = find(theNoisyPhi.begin(), theNoisyPhi.end(), iphi);
328  if (phi_iter != theNoisyPhi.end()) {
329  isAccepted = true;
330  }
331  return isAccepted;
332 }
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
std::vector< int > theNoisyPhi
string HPDNoiseLibraryReader::itos ( int  i)
private

Definition at line 385 of file HPDNoiseLibraryReader.cc.

References i, and alignCSCRings::s.

Referenced by getNoiseData().

385  {
386  stringstream s;
387 
388  s << i;
389  return s.str();
390 }
int i
Definition: DBlmapReader.cc:9
void HPDNoiseLibraryReader::Rannor ( double &  a,
double &  b 
)
private

Definition at line 372 of file HPDNoiseLibraryReader.cc.

References alignCSCRings::r, theRandFlat, x, detailsBasic3DVector::y, and detailsBasic3DVector::z.

Referenced by getIonFeedbackNoise().

372  {
373  double r,
374  x,
375  y,
376  z;
377 
378  y = theRandFlat->shoot();
379  z = theRandFlat->shoot();
380  x = z * 6.28318530717958623;
381  r = TMath::Sqrt(-2 * TMath::Log(y));
382  a = r * TMath::Sin(x);
383  b = r * TMath::Cos(x);
384 }
float float float z
double b
Definition: hdecay.h:120
CLHEP::RandFlat * theRandFlat
double a
Definition: hdecay.h:121
Definition: DDAxes.h:10
void HPDNoiseLibraryReader::setRandomEngine ( )
protected

Definition at line 67 of file HPDNoiseLibraryReader.cc.

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

Referenced by HPDNoiseLibraryReader().

67  {
69  if (!rng.isAvailable()) {
70  throw cms::Exception("Configuration") << "HcalHPDNoiseLibrary requires the RandomNumberGeneratorService\n"
71  "which is not present in the configuration file. You must add the service\n"
72  "in the configuration file or remove the modules that require it.";
73  }
74  setRandomEngine(rng->getEngine());
75 }
bool isAvailable() const
Definition: Service.h:46
virtual CLHEP::HepRandomEngine & getEngine() const =0
Use this to get the random number engine, this is the only function most users should call...
void HPDNoiseLibraryReader::setRandomEngine ( CLHEP::HepRandomEngine &  engine)
protected

Definition at line 76 of file HPDNoiseLibraryReader.cc.

References theRandFlat, and theRandGaussQ.

76  {
77  if (theRandGaussQ)
78  delete theRandGaussQ;
79 
80  if (theRandFlat)
81  delete theRandFlat;
82 
83  theRandGaussQ = new CLHEP::RandGaussQ(engine);
84  theRandFlat = new CLHEP::RandFlat(engine);
85 }
CLHEP::RandGaussQ * theRandGaussQ
CLHEP::RandFlat * theRandFlat
void HPDNoiseLibraryReader::shuffleData ( int  timeSliceId,
float *&  data 
)
private

Definition at line 333 of file HPDNoiseLibraryReader.cc.

References createTree::dd, i, gen::k, and max().

Referenced by getBiasedNoisyHcalDetIds(), and getNoisyHcalDetIds().

334 {
335  if(timeSliceId == -1 || (timeSliceId>=10)) return;
336  //make a local copy of input data
337  float Data[10] = {0.,0.,0.,0.,0.,0.,0.,0.,0.,0.};
338  for(int i=0;i<10;++i){
339  Data[i] = data[i];
340  }
341  int ts_max = -1;
342  float max = -999.;
343  for(int i=0;i<10;++i){
344  if(Data[i]>max){
345  max = data[i];
346  ts_max = i;
347  }
348  }
349  if((ts_max == -1)){//couldn't find ts_max, return the same value.
350  return;
351  }else{
352  // always shift the noise to the right by putting zeroes to the previous slices.
353  // the noise is pedestal subtracted. 0 value is acceptable.
354  int k = -1;
355  for(int i=0;i<10;++i){
356  data[i] = 0.;
357  int newIdx = timeSliceId+k;
358  float dd = 0.;
359  if(newIdx < 10){
360  data[newIdx] = Data[ts_max+k];
361  dd = Data[ts_max+k];
362  i = newIdx;
363  }
364  data[i] = dd;
365  ++k;
366  }
367 
368  }
369 }
int i
Definition: DBlmapReader.cc:9
const T & max(const T &a, const T &b)
int k[5][pyjets_maxn]
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82

Member Data Documentation

std::vector<float> HPDNoiseLibraryReader::theDischargeNoiseRate

Definition at line 96 of file HPDNoiseLibraryReader.h.

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

std::string HPDNoiseLibraryReader::theHPDName

Definition at line 104 of file HPDNoiseLibraryReader.h.

Referenced by getNoiseData(), and HPDNoiseLibraryReader().

std::vector<float> HPDNoiseLibraryReader::theIonFeedbackFirstPeakRate

Definition at line 97 of file HPDNoiseLibraryReader.h.

Referenced by fillRates(), and getIonFeedbackNoise().

std::vector<float> HPDNoiseLibraryReader::theIonFeedbackSecondPeakRate

Definition at line 98 of file HPDNoiseLibraryReader.h.

Referenced by fillRates(), and getIonFeedbackNoise().

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

Definition at line 103 of file HPDNoiseLibraryReader.h.

Referenced by fillRates(), and HPDNoiseLibraryReader().

std::vector<int> HPDNoiseLibraryReader::theNoisyPhi
CLHEP::RandFlat* HPDNoiseLibraryReader::theRandFlat
CLHEP::RandGaussQ* HPDNoiseLibraryReader::theRandGaussQ
HPDNoiseReader* HPDNoiseLibraryReader::theReader

Definition at line 102 of file HPDNoiseLibraryReader.h.

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