Typedefs | Functions
WireCell::Waveform Namespace Reference

Typedefs

typedef float real_t
 The type for the signal in each bin. More...
 
typedef std::complex< float > complex_t
 The type for the spectrum in each bin. More...
 
template<typename Val >
using Sequence = std::vector< Val >
 
typedef Sequence< real_trealseq_t
 
typedef Sequence< complex_tcompseq_t
 A complex-valued sequence, eg for discrete spectrum powers. More...
 
typedef std::pair< int, int > BinRange
 A half-open range of bins (from first bin to one past last bin) More...
 
typedef std::vector< BinRangeBinRangeList
 A list of bin ranges. More...
 
typedef std::map< int, BinRangeListChannelMasks
 Map channel number to a vector of BinRanges. More...
 
typedef std::map< std::string, ChannelMasksChannelMaskMap
 Collect channel masks by some label. More...
 
typedef std::pair< double, double > Period
 A range of time. More...
 
typedef std::pair< double, double > Band
 A range of frequency. More...
 
typedef std::pair< double, double > Domain
 

Functions

BinRangeList merge (const BinRangeList &br)
 Return a new list with any overlaps formed into unions. More...
 
BinRangeList merge (const BinRangeList &br1, const BinRangeList &br2)
 Merge two bin range lists, forming a union from any overlapping ranges. More...
 
ChannelMasks merge (const ChannelMasks &one, const ChannelMasks &two)
 Return a new mapping which is the union of all same channel masks. More...
 
void merge (ChannelMaskMap &one, ChannelMaskMap &two, std::map< std::string, std::string > &name_map)
 
int sample_count (const Domain &domain, double width)
 Return the number of samples needed to cover the domain with sample size width. More...
 
double sample_width (const Domain &domain, int count)
 Return the sample size if domain is equally sampled with given number of samples. More...
 
std::pair< int, int > sub_sample (const Domain &domain, int nsamples, const Domain &subdomain)
 Return the begin/end sample numbers inside the a subdomain of a domain with nsamples total. More...
 
template<typename Val >
Sequence< Val > resample (const Sequence< Val > &wave, const Domain &domain, int nsamples, const Domain &newdomain)
 
realseq_t real (const compseq_t &seq)
 Return the real part of the sequence. More...
 
realseq_t imag (const compseq_t &seq)
 Return the imaginary part of the sequence. More...
 
realseq_t magnitude (const compseq_t &seq)
 Return the magnitude or absolute value of the sequence. More...
 
realseq_t phase (const compseq_t &seq)
 Return the phase or arg part of the sequence. More...
 
template<typename Val >
void increase (Sequence< Val > &seq, Val scalar)
 Increase (shift) sequence values by scalar. More...
 
void increase (Sequence< float > &seq, double scalar)
 
template<typename Val >
void increase (Sequence< Val > &seq, const Sequence< Val > &other)
 Increase (shift) sequence values by values in another sequence. More...
 
template<typename Val >
void scale (Sequence< Val > &seq, Val scalar)
 Scale (multiply) sequence values by scalar. More...
 
void scale (Sequence< float > &seq, double scalar)
 
template<typename Val >
void scale (Sequence< Val > &seq, const Sequence< Val > &other)
 Scale (multiply) seq values by values from the other sequence. More...
 
template<typename Val >
void shrink (Sequence< Val > &seq, const Sequence< Val > &other)
 Shrink (divide) seq values by values from the other sequence. More...
 
std::pair< int, int > edge (const realseq_t &wave)
 
template<typename Val >
Val sum (const Sequence< Val > &seq)
 Return sum of all entries in sequence. More...
 
template<typename Val >
Val sum2 (const Sequence< Val > &seq)
 Return sum of square of all entries in sequence. More...
 
std::pair< double, double > mean_rms (const realseq_t &wave)
 
real_t median (realseq_t &wave)
 
real_t median_binned (realseq_t &wave)
 
real_t percentile (realseq_t &wave, real_t percentage)
 
real_t percentile_binned (realseq_t &wave, real_t percentage)
 
compseq_t dft (realseq_t seq)
 
realseq_t linear_convolve (Waveform::realseq_t in1, Waveform::realseq_t in2, bool truncate=true)
 
realseq_t replace_convolve (Waveform::realseq_t wave, Waveform::realseq_t newres, Waveform::realseq_t oldres, bool truncate=true)
 
realseq_t idft (compseq_t spec)
 
short most_frequent (const std::vector< short > &vals)
 Return the smallest, most frequent value to appear in vector. More...
 

Typedef Documentation

typedef std::pair<double,double> WireCell::Waveform::Band

A range of frequency.

Definition at line 69 of file Waveform.h.

typedef std::pair<int,int> WireCell::Waveform::BinRange

A half-open range of bins (from first bin to one past last bin)

Definition at line 38 of file Waveform.h.

A list of bin ranges.

Definition at line 41 of file Waveform.h.

Collect channel masks by some label.

Definition at line 59 of file Waveform.h.

Map channel number to a vector of BinRanges.

Definition at line 50 of file Waveform.h.

typedef std::complex<float> WireCell::Waveform::complex_t

The type for the spectrum in each bin.

Definition at line 21 of file Waveform.h.

A complex-valued sequence, eg for discrete spectrum powers.

Definition at line 34 of file Waveform.h.

typedef std::pair<double,double> WireCell::Waveform::Domain

A domain of a sequence is bounded by the time or frequency at the start of its first element and that at the end of its last.

Definition at line 75 of file Waveform.h.

typedef std::pair<double,double> WireCell::Waveform::Period

A range of time.

Definition at line 66 of file Waveform.h.

The type for the signal in each bin.

Definition at line 18 of file Waveform.h.

Definition at line 31 of file Waveform.h.

template<typename Val >
using WireCell::Waveform::Sequence = typedef std::vector<Val>

A sequence is an ordered array of values (real or complex). By itself it is not associated with a domain.

Definition at line 27 of file Waveform.h.

Function Documentation

Waveform::compseq_t WireCell::Waveform::dft ( realseq_t  seq)

Discrete Fourier transform of real sequence. Returns full spectrum. No normalization scaling applied

Definition at line 141 of file Waveform.cxx.

142 {
143  auto v = Eigen::Map<Eigen::VectorXf>(wave.data(), wave.size());
144  Eigen::FFT<Waveform::real_t> trans;
145  Eigen::VectorXcf ret = trans.fwd(v);
146  return compseq_t(ret.data(), ret.data()+ret.size());
147 }
Sequence< complex_t > compseq_t
A complex-valued sequence, eg for discrete spectrum powers.
Definition: Waveform.h:34
std::pair< int, int > WireCell::Waveform::edge ( const realseq_t wave)

Return a pair of indices into wave which bound non-zero region. First index is of first non-zero sample, second index is one past last non-zero sample. If entire wave is empty then both are set to size().

Definition at line 121 of file Waveform.cxx.

122 {
123  const int size = wave.size();
124  int imin=size, imax=size;
125 
126  for (int ind=0; ind < size; ++ind) {
127  const real_t val = wave[ind];
128  if (val != 0.0) {
129  if (imin == size) { // found start edge
130  imin = ind;
131  }
132  if (imin < size) {
133  imax = ind+1;
134  }
135  }
136  }
137  return std::make_pair(imin, imax);
138 }
float real_t
The type for the signal in each bin.
Definition: Waveform.h:18
decltype(auto) constexpr size(T &&obj)
ADL-aware version of std::size.
Definition: StdUtils.h:87
Waveform::realseq_t WireCell::Waveform::idft ( compseq_t  spec)

Inverse, discrete Fourier transform. Expects full spectrum (twice Nyquist frequency). Applies the 1/Nsamples normalization.

Definition at line 149 of file Waveform.cxx.

150 {
151  Eigen::FFT<Waveform::real_t> trans;
152  auto v = Eigen::Map<Eigen::VectorXcf>(spec.data(), spec.size());
153  Eigen::VectorXf ret;
154  trans.inv(ret, v);
155  return realseq_t(ret.data(), ret.data()+ret.size());
156 }
Sequence< real_t > realseq_t
Definition: Waveform.h:31
Waveform::realseq_t WireCell::Waveform::imag ( const compseq_t seq)

Return the imaginary part of the sequence.

Definition at line 60 of file Waveform.cxx.

61 {
62  return c2r(seq, [](Waveform::complex_t c) { return std::imag(c); });
63 }
realseq_t imag(const compseq_t &seq)
Return the imaginary part of the sequence.
Definition: Waveform.cxx:60
std::complex< float > complex_t
The type for the spectrum in each bin.
Definition: Waveform.h:21
Waveform::realseq_t c2r(const Waveform::compseq_t &seq, Func func)
Definition: Waveform.cxx:48
template<typename Val >
void WireCell::Waveform::increase ( Sequence< Val > &  seq,
Val  scalar 
)

Increase (shift) sequence values by scalar.

Definition at line 129 of file Waveform.h.

129  {
130  std::transform(seq.begin(), seq.end(), seq.begin(),
131  [scalar](Val x) { return x+scalar; });
132  }
tagset_t transform(const tagset_t &ts, const ruleset_t &rs, bool all_rules=true)
Definition: TagRules.cxx:43
list x
Definition: train.py:276
void WireCell::Waveform::increase ( Sequence< float > &  seq,
double  scalar 
)
inline

Definition at line 133 of file Waveform.h.

133  {
134  increase(seq, (float)scalar);
135  }
void increase(Sequence< Val > &seq, const Sequence< Val > &other)
Increase (shift) sequence values by values in another sequence.
Definition: Waveform.h:139
template<typename Val >
void WireCell::Waveform::increase ( Sequence< Val > &  seq,
const Sequence< Val > &  other 
)

Increase (shift) sequence values by values in another sequence.

Definition at line 139 of file Waveform.h.

139  {
140  std::transform(seq.begin(), seq.end(), other.begin(), seq.begin(),
141  std::plus<Val>());
142  }
tagset_t transform(const tagset_t &ts, const ruleset_t &rs, bool all_rules=true)
Definition: TagRules.cxx:43
Waveform::realseq_t WireCell::Waveform::linear_convolve ( Waveform::realseq_t  in1,
Waveform::realseq_t  in2,
bool  truncate = true 
)

Definition at line 159 of file Waveform.cxx.

162 {
163  size_t n1_orig = in1.size(), n2_orig = in2.size();
164  size_t n_out = n1_orig + n2_orig - 1;
165 
166  in1.resize(n_out, 0);
167  in2.resize(n_out, 0);
168 
169  auto v1 = Eigen::Map<Eigen::VectorXf>(in1.data(), in1.size());
170  auto v2 = Eigen::Map<Eigen::VectorXf>(in2.data(), in2.size());
171 
172  Eigen::FFT<Waveform::real_t> trans;
173 
174  Eigen::VectorXcf s1 = trans.fwd(v1);
175  Eigen::VectorXcf s2 = trans.fwd(v2);
176  Eigen::VectorXcf s12 = (s1.array() * s2.array()).matrix();
177  Eigen::VectorXf vret;
178  trans.inv(vret, s12);
179  realseq_t ret(vret.data(), vret.data()+vret.size());
180  if (truncate) {
181  ret.resize(n1_orig);
182  }
183  return ret;
184 }
Sequence< real_t > realseq_t
Definition: Waveform.h:31
Waveform::realseq_t WireCell::Waveform::magnitude ( const compseq_t seq)

Return the magnitude or absolute value of the sequence.

Definition at line 65 of file Waveform.cxx.

66 {
67  return c2r(seq, [](Waveform::complex_t c) { return std::abs(c); });
68 }
T abs(T value)
std::complex< float > complex_t
The type for the spectrum in each bin.
Definition: Waveform.h:21
Waveform::realseq_t c2r(const Waveform::compseq_t &seq, Func func)
Definition: Waveform.cxx:48
std::pair< double, double > WireCell::Waveform::mean_rms ( const realseq_t wave)

Definition at line 24 of file Waveform.cxx.

25 {
26  const int n = wf.size();
27  if (n==0) {
28  return std::make_pair<double,double>(0,0);
29  }
30  if (n==1) {
31  return std::make_pair<double,double>(wf[0],0);
32  }
33 
34  // if left as float, numerical precision will lead to many NaN for
35  // the RMS due to sometimes subtracting similar sized numbers.
36  std::vector<double> wfd(wf.begin(), wf.end());
37 
38  const double wsum = Waveform::sum(wfd);
39  const double w2sum = Waveform::sum2(wfd);
40  const double mean = wsum/n;
41  const double rms = sqrt( (w2sum - wsum*wsum/n) / n );
42  // const double rms = sqrt( (w2sum - wsum*wsum/n) / (n-1) );
43  return std::make_pair(mean,rms);
44 }
double rms(sqlite3 *db, std::string const &table_name, std::string const &column_name)
Definition: statistics.cc:39
Val sum2(const Sequence< Val > &seq)
Return sum of square of all entries in sequence.
Definition: Waveform.h:184
IFrame::pointer sum(std::vector< IFrame::pointer > frames, int ident)
Definition: FrameUtil.cxx:15
std::size_t n
Definition: format.h:3399
double mean(sqlite3 *db, std::string const &table_name, std::string const &column_name)
Definition: statistics.cc:15
Waveform::real_t WireCell::Waveform::median ( Waveform::realseq_t wave)

Definition at line 76 of file Waveform.cxx.

77 {
78  return percentile(wave,0.5);
79 }
real_t percentile(realseq_t &wave, real_t percentage)
Definition: Waveform.cxx:87
Waveform::real_t WireCell::Waveform::median_binned ( Waveform::realseq_t wave)

Definition at line 81 of file Waveform.cxx.

82 {
83  return percentile_binned(wave,0.5);
84 }
real_t percentile_binned(realseq_t &wave, real_t percentage)
Definition: Waveform.cxx:93
WireCell::Waveform::BinRangeList WireCell::Waveform::merge ( const BinRangeList br)

Return a new list with any overlaps formed into unions.

Definition at line 222 of file Waveform.cxx.

223 {
224  WireCell::Waveform::BinRangeList tmp(brl.begin(), brl.end());
226  sort(tmp.begin(), tmp.end());
227  Waveform::BinRange last_br = tmp[0];
228  out.push_back(last_br);
229 
230  for (size_t ind=1; ind<tmp.size(); ++ind) {
231  Waveform::BinRange this_br = tmp[ind];
232  if (out.back().second >= this_br.first) {
233  out.back().second = this_br.second;
234  continue;
235  }
236  out.push_back(this_br);
237  }
238  return out;
239 }
std::pair< int, int > BinRange
A half-open range of bins (from first bin to one past last bin)
Definition: Waveform.h:38
std::vector< BinRange > BinRangeList
A list of bin ranges.
Definition: Waveform.h:41
string tmp
Definition: languages.py:63
WireCell::Waveform::BinRangeList WireCell::Waveform::merge ( const BinRangeList br1,
const BinRangeList br2 
)

Merge two bin range lists, forming a union from any overlapping ranges.

Definition at line 243 of file Waveform.cxx.

245 {
247  out.reserve(br1.size() + br2.size());
248  out.insert(out.end(), br1.begin(), br1.end());
249  out.insert(out.end(), br2.begin(), br2.end());
250  return merge(out);
251 }
BinRangeList merge(const BinRangeList &br)
Return a new list with any overlaps formed into unions.
Definition: Waveform.cxx:222
std::vector< BinRange > BinRangeList
A list of bin ranges.
Definition: Waveform.h:41
WireCell::Waveform::ChannelMasks WireCell::Waveform::merge ( const ChannelMasks one,
const ChannelMasks two 
)

Return a new mapping which is the union of all same channel masks.

Definition at line 258 of file Waveform.cxx.

260 {
262  for (auto const &it : two) {
263  int ch = it.first;
264  out[ch] = merge(out[ch], it.second);
265  }
266  return out;
267 }
std::map< int, BinRangeList > ChannelMasks
Map channel number to a vector of BinRanges.
Definition: Waveform.h:50
BinRangeList merge(const BinRangeList &br)
Return a new list with any overlaps formed into unions.
Definition: Waveform.cxx:222
void WireCell::Waveform::merge ( ChannelMaskMap one,
ChannelMaskMap two,
std::map< std::string, std::string > &  name_map 
)

Definition at line 270 of file Waveform.cxx.

270  {
271 
272  // loop over second map
273  for (auto const& it: two){
274  std::string name = it.first;
275  std::string mapped_name;
276  if (name_map.find(name)!=name_map.end()){
277  mapped_name = name_map[name];
278  }else{
279  mapped_name = name;
280  }
281  if (one.find(mapped_name) != one.end()){
282  one[mapped_name] = merge(one[mapped_name],it.second);
283  }else{
284  one[mapped_name] = it.second;
285  }
286  }
287 }
static QCString name
Definition: declinfo.cpp:673
std::string string
Definition: nybbler.cc:12
BinRangeList merge(const BinRangeList &br)
Return a new list with any overlaps formed into unions.
Definition: Waveform.cxx:222
short WireCell::Waveform::most_frequent ( const std::vector< short > &  vals)

Return the smallest, most frequent value to appear in vector.

Definition at line 289 of file Waveform.cxx.

290 {
291  const size_t nbins = 1<<16;
292  std::vector<unsigned int> hist(nbins, 0);
293  for (unsigned short val : vals) {
294  hist[val] += 1;
295  }
296  auto it = std::max_element(hist.begin(), hist.end());
297  return it - hist.begin(); // casts back to signed short
298 }
Waveform::real_t WireCell::Waveform::percentile ( Waveform::realseq_t wave,
real_t  percentage 
)

Definition at line 87 of file Waveform.cxx.

88 {
89  std::nth_element(wave.begin(), wave.begin()+wave.size()*percentage, wave.end());
90  return wave.at(wave.size()*percentage);
91 }
Waveform::real_t WireCell::Waveform::percentile_binned ( Waveform::realseq_t wave,
real_t  percentage 
)

Definition at line 93 of file Waveform.cxx.

93  {
94  const auto mm = std::minmax_element(wave.begin(), wave.end());
95  const auto vmin = *mm.first;
96  const auto vmax = *mm.second;
97  const int nbins = wave.size();
98  const auto binsize = (vmax-vmin)/nbins;
100  for (auto val : wave) {
101  int bin = int(round((val - vmin)/binsize));
102  bin = std::max(0, bin);
103  bin = std::min(nbins-1, bin);
104  hist[bin] += 1.0;
105  }
106 
107  const int imed = wave.size() * percentage;
108  int count = 0;
109  for (int ind=0; ind<nbins; ++ind) {
110  count += hist[ind];
111  if (count > imed) {
112  float ret = vmin + ind*binsize;
113  return ret;
114  }
115  }
116  // can't reach here, return bogus value.
117  return vmin + (vmax-vmin)*percentage;
118 }
Sequence< real_t > realseq_t
Definition: Waveform.h:31
static int max(int a, int b)
static const double mm
Definition: Units.h:73
T min(sqlite3 *const db, std::string const &table_name, std::string const &column_name)
Definition: statistics.h:55
QTextStream & bin(QTextStream &s)
Waveform::realseq_t WireCell::Waveform::phase ( const compseq_t seq)

Return the phase or arg part of the sequence.

Definition at line 70 of file Waveform.cxx.

71 {
72  return c2r(seq, [](Waveform::complex_t c) { return std::arg(c); });
73 }
internal::named_arg< T, char > arg(string_view name, const T &arg)
Definition: core.h:1391
std::complex< float > complex_t
The type for the spectrum in each bin.
Definition: Waveform.h:21
Waveform::realseq_t c2r(const Waveform::compseq_t &seq, Func func)
Definition: Waveform.cxx:48
Waveform::realseq_t WireCell::Waveform::real ( const compseq_t seq)

Return the real part of the sequence.

Definition at line 55 of file Waveform.cxx.

56 {
57  return c2r(seq, [](Waveform::complex_t c) { return std::real(c); });
58 }
realseq_t real(const compseq_t &seq)
Return the real part of the sequence.
Definition: Waveform.cxx:55
std::complex< float > complex_t
The type for the spectrum in each bin.
Definition: Waveform.h:21
Waveform::realseq_t c2r(const Waveform::compseq_t &seq, Func func)
Definition: Waveform.cxx:48
Waveform::realseq_t WireCell::Waveform::replace_convolve ( Waveform::realseq_t  wave,
Waveform::realseq_t  newres,
Waveform::realseq_t  oldres,
bool  truncate = true 
)

Definition at line 187 of file Waveform.cxx.

191 {
192  size_t sizes[3] = {wave.size(), newres.size(), oldres.size()};
193  size_t n_out = sizes[0]+sizes[1]+sizes[2] - *std::min_element(sizes, sizes+3) - 1;
194 
195  wave.resize(n_out, 0);
196  newres.resize(n_out, 0);
197  oldres.resize(n_out, 0);
198 
199  auto v1 = Eigen::Map<Eigen::VectorXf>(wave.data(), wave.size());
200  auto v2 = Eigen::Map<Eigen::VectorXf>(newres.data(), newres.size());
201  auto v3 = Eigen::Map<Eigen::VectorXf>(oldres.data(), oldres.size());
202 
203  Eigen::FFT<Waveform::real_t> trans;
204 
205  Eigen::VectorXcf s1 = trans.fwd(v1);
206  Eigen::VectorXcf s2 = trans.fwd(v2);
207  Eigen::VectorXcf s3 = trans.fwd(v3);
208 
209  Eigen::VectorXcf s123 = (s1.array() * s2.array() / s3.array()).matrix();
210 
211  Eigen::VectorXf vret;
212  trans.inv(vret, s123);
213  realseq_t ret(vret.data(), vret.data()+vret.size());
214  if (truncate) {
215  ret.resize(sizes[0]);
216  }
217  return ret;
218 }
Sequence< real_t > realseq_t
Definition: Waveform.h:31
template<typename Val >
Sequence<Val> WireCell::Waveform::resample ( const Sequence< Val > &  wave,
const Domain domain,
int  nsamples,
const Domain newdomain 
)

Return a new sequence resampled and interpolated from the original wave defined over the domain to a new domain of nsamples.

Definition at line 92 of file Waveform.h.

92  {
93  const int oldnsamples = wave.size();
94  const double oldstep = sample_width(domain, oldnsamples);
95  const double step = sample_width(newdomain, nsamples);
96  Sequence<Val> ret;
97  for (int ind=0; ind<nsamples; ++ind) {
98  double cursor = newdomain.first + ind*step;
99  double oldfracsteps = (cursor-domain.first)/oldstep;
100  int oldind = int(oldfracsteps);
101  if (cursor <= domain.first || oldind <= 0) {
102  ret.push_back(wave[0]);
103  continue;
104  }
105  if (cursor >= domain.second or oldind+1 >= oldnsamples) {
106  ret.push_back(wave[oldnsamples-1]);
107  continue;
108  }
109  double d1 = oldfracsteps - oldstep*oldind;
110  double d2 = oldstep - d1;
111  Val newval = (wave[oldind] * d1 + wave[oldind+1]*d2) / oldstep;
112  ret.push_back(newval);
113  }
114  return ret;
115  }
double sample_width(const Domain &domain, int count)
Return the sample size if domain is equally sampled with given number of samples. ...
Definition: Waveform.h:81
int WireCell::Waveform::sample_count ( const Domain domain,
double  width 
)
inline

Return the number of samples needed to cover the domain with sample size width.

Definition at line 77 of file Waveform.h.

77  {
78  return (domain.second-domain.first)/width;
79  }
const double width
double WireCell::Waveform::sample_width ( const Domain domain,
int  count 
)
inline

Return the sample size if domain is equally sampled with given number of samples.

Definition at line 81 of file Waveform.h.

81  {
82  return (domain.second-domain.first)/count;
83  }
template<typename Val >
void WireCell::Waveform::scale ( Sequence< Val > &  seq,
Val  scalar 
)

Scale (multiply) sequence values by scalar.

Definition at line 146 of file Waveform.h.

146  {
147  std::transform(seq.begin(), seq.end(), seq.begin(),
148  [scalar](Val x) { return x*scalar; });
149  }
tagset_t transform(const tagset_t &ts, const ruleset_t &rs, bool all_rules=true)
Definition: TagRules.cxx:43
list x
Definition: train.py:276
void WireCell::Waveform::scale ( Sequence< float > &  seq,
double  scalar 
)
inline

Definition at line 150 of file Waveform.h.

150  {
151  scale(seq, (float)scalar);
152  }
void scale(Sequence< Val > &seq, const Sequence< Val > &other)
Scale (multiply) seq values by values from the other sequence.
Definition: Waveform.h:156
template<typename Val >
void WireCell::Waveform::scale ( Sequence< Val > &  seq,
const Sequence< Val > &  other 
)

Scale (multiply) seq values by values from the other sequence.

Definition at line 156 of file Waveform.h.

156  {
157  std::transform(seq.begin(), seq.end(), other.begin(), seq.begin(),
158  std::multiplies<Val>());
159  }
tagset_t transform(const tagset_t &ts, const ruleset_t &rs, bool all_rules=true)
Definition: TagRules.cxx:43
template<typename Val >
void WireCell::Waveform::shrink ( Sequence< Val > &  seq,
const Sequence< Val > &  other 
)

Shrink (divide) seq values by values from the other sequence.

Definition at line 162 of file Waveform.h.

162  {
163  std::transform(seq.begin(), seq.end(), other.begin(), seq.begin(),
164  std::divides<Val>());
165  }
tagset_t transform(const tagset_t &ts, const ruleset_t &rs, bool all_rules=true)
Definition: TagRules.cxx:43
std::pair< int, int > WireCell::Waveform::sub_sample ( const Domain domain,
int  nsamples,
const Domain subdomain 
)

Return the begin/end sample numbers inside the a subdomain of a domain with nsamples total.

Definition at line 15 of file Waveform.cxx.

16 {
17  const double bin = sample_width(domain, nsamples);
18  int beg = ceil((subdomain.first - domain.first) / bin);
19  int end = nsamples- ceil((domain.second - subdomain.second) / bin);
20  return std::make_pair(std::max(0,beg), std::min(nsamples, end));
21 }
decltype(auto) constexpr end(T &&obj)
ADL-aware version of std::end.
Definition: StdUtils.h:72
double sample_width(const Domain &domain, int count)
Return the sample size if domain is equally sampled with given number of samples. ...
Definition: Waveform.h:81
static int max(int a, int b)
T min(sqlite3 *const db, std::string const &table_name, std::string const &column_name)
Definition: statistics.h:55
QTextStream & bin(QTextStream &s)
template<typename Val >
Val WireCell::Waveform::sum ( const Sequence< Val > &  seq)

Return sum of all entries in sequence.

Definition at line 178 of file Waveform.h.

178  {
179  return std::accumulate(seq.begin(), seq.end(), Val());
180  }
template<typename Val >
Val WireCell::Waveform::sum2 ( const Sequence< Val > &  seq)

Return sum of square of all entries in sequence.

Definition at line 184 of file Waveform.h.

184  {
185  return std::accumulate(seq.begin(), seq.end(), Val(),
186  [](const Val& bucket, Val x) {
187  return bucket + x*x;
188  });
189  }
list x
Definition: train.py:276