CMS 3D CMS Logo

List of all members | Public Member Functions | Protected Member Functions | Private Attributes
GaussianTailNoiseGenerator Class Reference

#include <GaussianTailNoiseGenerator.h>

Public Member Functions

 GaussianTailNoiseGenerator ()
 
void generate (int NumberOfchannels, float threshold, float noiseRMS, std::map< int, float > &theMap, CLHEP::HepRandomEngine *)
 
void generate (int NumberOfchannels, float threshold, float noiseRMS, std::vector< std::pair< int, float >> &, CLHEP::HepRandomEngine *)
 
void generateRaw (float noiseRMS, std::vector< double > &, CLHEP::HepRandomEngine *)
 

Protected Member Functions

double generate_gaussian_tail (const double, const double, CLHEP::HepRandomEngine *)
 
int * getRandomChannels (int, int, CLHEP::HepRandomEngine *)
 

Private Attributes

int channel512_ [512]
 
int channel768_ [768]
 

Detailed Description

Generation of random noise on "numberOfChannels" channels with a given threshold. The generated noise :

Initial author : Veronique Lefebure 08.10.98
according to the FORTRAN code tgreset.F from Pascal Vanlaer
Modified by C. Delaere 01.10.09

Fills in a map < channel number, generated noise >

Definition at line 24 of file GaussianTailNoiseGenerator.h.

Constructor & Destructor Documentation

◆ GaussianTailNoiseGenerator()

GaussianTailNoiseGenerator::GaussianTailNoiseGenerator ( )

Definition at line 11 of file GaussianTailNoiseGenerator.cc.

References channel512_, channel768_, and mps_fire::i.

11  {
12  // we have two cases: 512 and 768 channels
13  // other cases are not allowed so far (performances issue)
14  for (unsigned int i = 0; i < 512; ++i)
15  channel512_[i] = i;
16  for (unsigned int i = 0; i < 768; ++i)
17  channel768_[i] = i;
18 }

Member Function Documentation

◆ generate() [1/2]

void GaussianTailNoiseGenerator::generate ( int  NumberOfchannels,
float  threshold,
float  noiseRMS,
std::map< int, float > &  theMap,
CLHEP::HepRandomEngine *   
)

◆ generate() [2/2]

void GaussianTailNoiseGenerator::generate ( int  NumberOfchannels,
float  threshold,
float  noiseRMS,
std::vector< std::pair< int, float >> &  theVector,
CLHEP::HepRandomEngine *  engine 
)

Definition at line 55 of file GaussianTailNoiseGenerator.cc.

References DMR_cfg::cerr, ewkTauDQM_cfi::channels, generate_gaussian_tail(), getRandomChannels(), mps_fire::i, hgchebackDigitizer_cfi::noise, mps_fire::result, mps_update::status, and remoteMonitoring_LASER_era2018_cfg::threshold.

59  {
60  // Compute number of channels with noise above threshold
61  // Gaussian tail probability
62  gsl_sf_result result;
63  int status = gsl_sf_erf_Q_e(threshold, &result);
64  if (status != 0)
65  std::cerr << "GaussianTailNoiseGenerator::could not compute gaussian tail "
66  "probability for the threshold chosen"
67  << std::endl;
68  double probabilityLeft = result.val;
69  double meanNumberOfNoisyChannels = probabilityLeft * NumberOfchannels;
70 
71  CLHEP::RandPoissonQ randPoissonQ(*engine, meanNumberOfNoisyChannels);
72  int numberOfNoisyChannels = randPoissonQ.fire();
73 
74  if (numberOfNoisyChannels > NumberOfchannels)
75  numberOfNoisyChannels = NumberOfchannels;
76 
77  // Compute the list of noisy channels
78  theVector.reserve(numberOfNoisyChannels);
79  float lowLimit = threshold * noiseRMS;
80  int *channels = getRandomChannels(numberOfNoisyChannels, NumberOfchannels, engine);
81 
82  for (int i = 0; i < numberOfNoisyChannels; i++) {
83  // Find random noise value
84  double noise = generate_gaussian_tail(lowLimit, noiseRMS, engine);
85  // Fill in the vector
86  theVector.push_back(std::pair<int, float>(channels[i], noise));
87  }
88 }
double generate_gaussian_tail(const double, const double, CLHEP::HepRandomEngine *)
int * getRandomChannels(int, int, CLHEP::HepRandomEngine *)

◆ generate_gaussian_tail()

double GaussianTailNoiseGenerator::generate_gaussian_tail ( const double  a,
const double  sigma,
CLHEP::HepRandomEngine *  engine 
)
protected

Definition at line 141 of file GaussianTailNoiseGenerator.cc.

References a, dqm-mbProfile::log, alignCSCRings::s, mathSSE::sqrt(), findQualityFiles::v, and x.

Referenced by generate().

143  {
144  /* Returns a gaussian random variable larger than a
145  * This implementation does one-sided upper-tailed deviates.
146  */
147 
148  double s = a / sigma;
149 
150  if (s < 1) {
151  /*
152  For small s, use a direct rejection method. The limit s < 1
153  can be adjusted to optimise the overall efficiency
154  */
155  double x;
156 
157  do {
158  x = CLHEP::RandGaussQ::shoot(engine, 0., 1.0);
159  } while (x < s);
160  return x * sigma;
161 
162  } else {
163  /* Use the "supertail" deviates from the last two steps
164  * of Marsaglia's rectangle-wedge-tail method, as described
165  * in Knuth, v2, 3rd ed, pp 123-128. (See also exercise 11, p139,
166  * and the solution, p586.)
167  */
168 
169  double u, v, x;
170 
171  do {
172  u = CLHEP::RandFlat::shoot(engine);
173  do {
174  v = CLHEP::RandFlat::shoot(engine);
175  } while (v == 0.0);
176  x = sqrt(s * s - 2 * log(v));
177  } while (x * u > s);
178  return x * sigma;
179  }
180 }
T sqrt(T t)
Definition: SSEVec.h:19
double a
Definition: hdecay.h:121

◆ generateRaw()

void GaussianTailNoiseGenerator::generateRaw ( float  noiseRMS,
std::vector< double > &  theVector,
CLHEP::HepRandomEngine *  engine 
)

Definition at line 107 of file GaussianTailNoiseGenerator.cc.

References mps_fire::i.

109  {
110  // it was shown that a complex approach, inspired from the ZS case,
111  // does not allow to gain much.
112  // A cut at 2 sigmas only saves 25% of the processing time, while the cost
113  // in terms of meaning is huge.
114  // We therefore use here the trivial approach (as in the early 2XX cycle)
115  unsigned int numberOfchannels = theVector.size();
116  for (unsigned int i = 0; i < numberOfchannels; ++i) {
117  if (theVector[i] == 0)
118  theVector[i] = CLHEP::RandGaussQ::shoot(engine, 0., noiseRMS);
119  }
120 }

◆ getRandomChannels()

int * GaussianTailNoiseGenerator::getRandomChannels ( int  numberOfNoisyChannels,
int  numberOfchannels,
CLHEP::HepRandomEngine *  engine 
)
protected

Definition at line 122 of file GaussianTailNoiseGenerator.cc.

References mps_check::array, b, channel512_, channel768_, createfilelist::int, and dqmiolumiharvest::j.

Referenced by generate().

124  {
125  if (numberOfNoisyChannels > numberOfchannels)
126  numberOfNoisyChannels = numberOfchannels;
127  int *array = channel512_;
128  if (numberOfchannels == 768)
129  array = channel768_;
130  int theChannelNumber;
131  for (int j = 0; j < numberOfNoisyChannels; ++j) {
132  theChannelNumber = (int)CLHEP::RandFlat::shoot(engine, numberOfchannels - j) + j;
133  // swap the two array elements... this is optimized by the compiler
134  int b = array[j];
135  array[j] = array[theChannelNumber];
136  array[theChannelNumber] = b;
137  }
138  return array;
139 }
double b
Definition: hdecay.h:120

Member Data Documentation

◆ channel512_

int GaussianTailNoiseGenerator::channel512_[512]
private

Definition at line 53 of file GaussianTailNoiseGenerator.h.

Referenced by GaussianTailNoiseGenerator(), and getRandomChannels().

◆ channel768_

int GaussianTailNoiseGenerator::channel768_[768]
private

Definition at line 54 of file GaussianTailNoiseGenerator.h.

Referenced by GaussianTailNoiseGenerator(), and getRandomChannels().