CMS 3D CMS Logo

imath.h
Go to the documentation of this file.
1 // // Integer representation of floating point arithmetic suitable for FPGA designs
2 //
3 // Author: Yuri Gershtein
4 // Date: March 2018
5 //
6 // Functionality:
7 //
8 // *note* all integers are assumed to be signed
9 //
10 // all variables have units, stored in a map <string,int>, with string a unit (i.e. "phi") and int the power
11 // "2" is always present in the map, and it's int pair is referred to as 'shift'
12 // units are properly combined / propagated through calculations
13 // adding/subtracting variables with different units throws an exception
14 // adding/subtracting variables with different shifts is allowed and is handled correctly
15 //
16 // calculate() method re-calculates the variable double and int values based on its operands
17 // returns false in case of overflows and/or mismatches between double and int calculations.
18 //
19 // the maximum and minimum values that the variable assumes are stored and updated each time calculate() is called
20 // if IMATH_ROOT is defined, all values are also stored in a histogram
21 //
22 // VarDef (string name, string units, double fmax, double K):
23 // define variable with bit value fval = K*ival, and maximum absolute value fmax.
24 // calculates nbins on its own
25 // one can assign value to it using set_ methods
26 //
27 // VarParam (string name, string units, double fval, int nbits):
28 // define a parameter. K is calculated based on the fval and nbits
29 //
30 // or (string name, std::string units, double fval, double K):
31 // define a parameer with bit value fval = K*ival.
32 // calculates nbins on its own
33 //
34 // VarAdd (string name, VarBase *p1, VarBase *p2, double range = -1, int nmax = 18):
35 // VarSubtract(string name, VarBase *p1, VarBase *p2, double range = -1, int nmax = 18):
36 // add/subtract variables. Bit length increases by 1, but capped at nmax
37 // if range>0 specified, bit length is decreased to drop unnecessary high bits
38 //
39 // VarMult (string name, VarBase *p1, VarBase *p2, double range = -1, int nmax = 18):
40 // multiplication. Bit length is a sum of the lengths of the operads, but capped at nmax
41 // if range>0 specified, bit length is decreased to drop unnecessary high bits or post-shift is reduced
42 //
43 // VarTimesC (string name, VarBase *p1, double cF, int ps = 17):
44 // multiplication by a constant. Bit length stays the same
45 // ps defines number of bits used to represent the constant
46 //
47 // VarDSPPostadd (string name, VarBase *p1, VarBase *p2, VarBase *p3, double range = -1, int nmax = 18):
48 // explicit instantiation of the 3-clock DSP postaddition: p1*p2+p3
49 // range and nmax have the same meaning as for the VarMult.
50 //
51 // VarShift (string name, VarBase *p1, int shift):
52 // shifts the variable right by shift (equivalent to multiplication by pow(2, -shift));
53 // Units stay the same, nbits are adjusted.
54 //
55 // VarShiftround (string name, VarBase *p1, int shift):
56 // shifts the variable right by shift, but doing rounding, i.e.
57 // (p>>(shift-1)+1)>>1;
58 // Units stay the same, nbits are adjusted.
59 //
60 // VarNeg (string name, VarBase *p1):
61 // multiplies the variable by -1
62 //
63 // VarInv (string name, VarBase *p1, double offset, int nbits, int n, unsigned int shift, mode m, int nbaddr=-1):
64 // LUT-based inversion, f = 1./(offset + f1) and i = 2^n / (offsetI + i1)
65 // nbits is the width of the LUT (signed)
66 // m is from enum mode {pos, neg, both} and refers to possible sign values of f
67 // for pos and neg, the most significant bit of p1 (i.e. the sign bit) is ignored
68 // shift is a shift applied in i1<->address conversions (used to reduce size of LUT)
69 // nbaddr: if not specified, it is taken to be equal to p1->nbits()
70 //
71 //
72 // VarNounits (string name, VarBase *p1, int ps = 17):
73 // convert a number with units to a number - needed for trig function expansion (i.e. 1 - 0.5*phi^2)
74 // ps is a number of bits to represent the unit conversion constant
75 //
76 // VarAdjustK (string name, VarBase *p1, double Knew, double epsilon = 1e-5, bool do_assert = false, int nbits = -1)
77 // adjust variable shift so the K is as close to Knew as possible (needed for bit length adjustments)
78 // if do_assert is true, throw an exeption if Knew/Kold is not a power of two
79 // epsilon is a comparison precision, nbits forces the bit length (possibly discarding MSBs)
80 //
81 // VarAdjustKR (string name, VarBase *p1, double Knew, double epsilon = 1e-5, bool do_assert = false, int nbits = -1)
82 // - same as adjustK(), but with rounding, and therefore latency=1
83 //
84 // bool calculate(int debug_level) runs through the entire formula tree recalculating both ineteger and floating point values
85 // returns true if everything is OK, false if obvious problems with the calculation exist, i.e
86 // - integer value does not fit into the alotted number of bins
87 // - integer value is more then 10% or more then 2 away from fval_/K_
88 // debug_level: 0 - no warnings
89 // 1 - limited warning
90 // 2 - as 1, but also include explicit warnings when LUT was used out of its range
91 // 3 - maximum complaints level
92 //
93 // VarFlag (string name, VarBase *cut_var, VarBase *...)
94 //
95 // flag to apply cuts defined for any variable. When output as Verilog, the flag
96 // is true if and only if the following conditions are all true:
97 // 1) the cut defined by each VarCut pointer in the argument list must be passed
98 // by the associated variable
99 // 2) each VarBase pointer in the argument list that is not also a VarCut
100 // pointer must pass all of its associated cuts
101 // 3) all children of the variables in the argument list must pass all of their
102 // associated cuts
103 // The VarFlag::passes() method replicates the behavior of the output Verilog,
104 // returning true if and only if the above conditions are all true. The
105 // VarBase::local_passes() method can be used to query if a given variable passes
106 // its associated cuts, regardless of whether its children do.
107 //
108 #ifndef L1Trigger_TrackFindingTracklet_interface_imath_h
109 #define L1Trigger_TrackFindingTracklet_interface_imath_h
110 
111 //use root if uncommented
112 //#ifndef CMSSW_GIT_HASH
113 //#define IMATH_ROOT
114 //#endif
115 
116 #include <limits>
117 #include <iostream>
118 #include <fstream>
119 #include <vector>
120 #include <map>
121 #include <cmath>
122 #include <sstream>
123 #include <string>
124 #include <cassert>
125 #include <set>
126 
130 
131 #ifdef IMATH_ROOT
132 #include <TH2F.h>
133 #include <TFile.h>
134 #include <TCanvas.h>
135 #include <TTree.h>
136 #endif
137 
138 //operation latencies for proper HDL pipelining
139 #define MULT_LATENCY 1
140 #define LUT_LATENCY 2
141 #define DSP_LATENCY 3
142 
143 // Print out information on the pass/fail status of all variables. Warning:
144 // this outputs a lot of information for busy events!
145 
146 namespace trklet {
147 
148  struct imathGlobals {
149  bool printCutInfo_{false};
150 #ifdef IMATH_ROOT
151  TFile *h_file_ = new TFile("imath.root", "RECREATE");
152  bool use_root;
153 #endif
154  };
155 
156  class VarCut;
157  class VarFlag;
158 
159  class VarBase {
160  public:
162  globals_ = globals;
163  p1_ = p1;
164  p2_ = p2;
165  p3_ = p3;
166  name_ = name;
167  latency_ = l;
168  int step1 = (p1) ? p1->step() + p1->latency() : 0;
169  int step2 = (p2) ? p2->step() + p2->latency() : 0;
170  step_ = std::max(step1, step2);
171 
172  cuts_.clear();
173  cut_var_ = nullptr;
174 
175  pipe_counter_ = 0;
176  pipe_delays_.clear();
177 
180  readytoprint_ = true;
181  readytoanalyze_ = true;
182  usedasinput_ = false;
183  Kmap_.clear();
184  Kmap_["2"] = 0; // initially, zero shift
185 #ifdef IMATH_ROOT
186  h_ = 0;
187  h_nbins_ = 1024;
188  h_precision_ = 0.02;
189 #endif
190  }
191  virtual ~VarBase() {
192 #ifdef IMATH_ROOT
193  if (globals_->h_file_) {
194  globals_->h_file_->ls();
195  globals_->h_file_->Close();
196  globals_->h_file_ = 0;
197  }
198 #endif
199  }
200 
201  static struct Verilog {
202  } verilog;
203  static struct HLS {
204  } hls;
205 
206  std::string kstring() const;
207  std::string name() const { return name_; }
208  std::string op() const { return op_; }
209  VarBase *p1() const { return p1_; }
210  VarBase *p2() const { return p2_; }
211  VarBase *p3() const { return p3_; }
212  double fval() const { return fval_; }
213  long int ival() const { return ival_; }
214 
215  bool local_passes() const;
216  void passes(std::map<const VarBase *, std::vector<bool> > &passes,
217  const std::map<const VarBase *, std::vector<bool> > *const previous_passes = nullptr) const;
218  void print_cuts(std::map<const VarBase *, std::set<std::string> > &cut_strings,
219  const int step,
220  Verilog,
221  const std::map<const VarBase *, std::set<std::string> > *const previous_cut_strings = nullptr) const;
222  void print_cuts(std::map<const VarBase *, std::set<std::string> > &cut_strings,
223  const int step,
224  HLS,
225  const std::map<const VarBase *, std::set<std::string> > *const previous_cut_strings = nullptr) const;
226  void add_cut(VarCut *cut, const bool call_set_cut_var = true);
227  VarBase *cut_var();
228 
229  double minval() const { return minval_; }
230  double maxval() const { return maxval_; }
231  void analyze();
232 #ifdef IMATH_ROOT
233  TH2F *h() { return h_; }
234 #endif
235  void reset() {
238 #ifdef IMATH_ROOT
239  h_->Clear();
240 #endif
241  }
242 
243  int nbits() const { return nbits_; }
244  std::map<std::string, int> Kmap() const { return Kmap_; }
245  double range() const { return (1 << (nbits_ - 1)) * K_; } // everything is signed
246  double K() const { return K_; };
247  int shift() const { return Kmap_.at("2"); }
248 
249  void makeready();
250  int step() const { return step_; }
251  int latency() const { return latency_; }
252  void add_latency(unsigned int l) { latency_ += l; } //only call before using the variable in calculation!
253  bool calculate(int debug_level);
254  bool calculate() { return calculate(0); }
255  virtual void local_calculate() {}
256  virtual void print(std::ofstream &fs, Verilog, int l1 = 0, int l2 = 0, int l3 = 0) {
257  fs << "// VarBase here. Soemthing is wrong!! " << l1 << ", " << l2 << ", " << l3 << "\n";
258  }
259  virtual void print(std::ofstream &fs, HLS, int l1 = 0, int l2 = 0, int l3 = 0) {
260  fs << "// VarBase here. Soemthing is wrong!! " << l1 << ", " << l2 << ", " << l3 << "\n";
261  }
262  void print_step(int step, std::ofstream &fs, Verilog);
263  void print_step(int step, std::ofstream &fs, HLS);
264  void print_all(std::ofstream &fs, Verilog);
265  void print_all(std::ofstream &fs, HLS);
266  void print_truncation(std::string &t, const std::string &o1, const int ps, Verilog) const;
267  void print_truncation(std::string &t, const std::string &o1, const int ps, HLS) const;
268  void inputs(std::vector<VarBase *> *vd); //collect all inputs
269 
270  int pipe_counter() { return pipe_counter_; }
272  void add_delay(int i) { pipe_delays_.push_back(i); }
273  bool has_delay(int i); //returns true if already have this variable delayed.
274  static void verilog_print(const std::vector<VarBase *> &v, std::ofstream &fs) { design_print(v, fs, verilog); }
275  static void hls_print(const std::vector<VarBase *> &v, std::ofstream &fs) { design_print(v, fs, hls); }
276  static void design_print(const std::vector<VarBase *> &v, std::ofstream &fs, Verilog);
277  static void design_print(const std::vector<VarBase *> &v, std::ofstream &fs, HLS);
278  static std::string pipe_delay(VarBase *v, int nbits, int delay);
279  std::string pipe_delays(const int step);
280  static std::string pipe_delay_wire(VarBase *v, std::string name_delayed, int nbits, int delay);
281 
282 #ifdef IMATH_ROOT
283  static TTree *addToTree(imathGlobals *globals, VarBase *v, char *s = 0);
284  static TTree *addToTree(imathGlobals *globals, int *v, char *s);
285  static TTree *addToTree(imathGlobals *globals, double *v, char *s);
286  static void fillTree(imathGlobals *globals);
287  static void writeTree(imathGlobals *globals);
288 #endif
289 
290  void dump_msg();
291  std::string dump();
292  static std::string itos(int i);
293 
294  protected:
300  std::string op_; // operation
301  int latency_; // number of clock cycles for the operation (for HDL output)
302  int step_; // step number in the calculation (for HDL output)
303 
304  double fval_; // exact calculation
305  long int ival_; // integer calculation
306  double val_; // integer calculation converted to double, ival_*K
307 
308  std::vector<VarBase *> cuts_;
310 
311  int nbits_;
312  double K_;
313  std::map<std::string, int> Kmap_;
314 
316  std::vector<int> pipe_delays_;
317 
321 
322  double minval_;
323  double maxval_;
324 #ifdef IMATH_ROOT
325  void set_hist_pars(int n = 256, double p = 0.05) {
326  h_nbins_ = n;
327  h_precision_ = p;
328  }
329  int h_nbins_;
330  double h_precision_;
331  TH2F *h_;
332 #endif
333  };
334 
335  class VarAdjustK : public VarBase {
336  public:
339  VarBase *p1,
340  double Knew,
341  double epsilon = 1e-5,
342  bool do_assert = false,
343  int nbits = -1)
344  : VarBase(globals, name, p1, nullptr, nullptr, 0) {
345  op_ = "adjustK";
346  K_ = p1->K();
347  Kmap_ = p1->Kmap();
348 
349  double r = Knew / K_;
350 
351  lr_ = (r > 1) ? log2(r) + epsilon : log2(r);
352  K_ = K_ * pow(2, lr_);
353  if (do_assert)
354  assert(std::abs(Knew / K_ - 1) < epsilon);
355 
356  if (nbits > 0)
357  nbits_ = nbits;
358  else
359  nbits_ = p1->nbits() - lr_;
360 
361  Kmap_["2"] = Kmap_["2"] + lr_;
362  }
363 
364  ~VarAdjustK() override = default;
365 
366  void adjust(double Knew, double epsilon = 1e-5, bool do_assert = false, int nbits = -1);
367 
368  void local_calculate() override;
369  void print(std::ofstream &fs, Verilog, int l1 = 0, int l2 = 0, int l3 = 0) override;
370  void print(std::ofstream &fs, HLS, int l1 = 0, int l2 = 0, int l3 = 0) override;
371 
372  protected:
373  int lr_;
374  };
375 
376  class VarAdjustKR : public VarBase {
377  public:
380  VarBase *p1,
381  double Knew,
382  double epsilon = 1e-5,
383  bool do_assert = false,
384  int nbits = -1)
385  : VarBase(globals, name, p1, nullptr, nullptr, 1) {
386  op_ = "adjustKR";
387  K_ = p1->K();
388  Kmap_ = p1->Kmap();
389 
390  double r = Knew / K_;
391 
392  lr_ = (r > 1) ? log2(r) + epsilon : log2(r);
393  K_ = K_ * pow(2, lr_);
394  if (do_assert)
395  assert(std::abs(Knew / K_ - 1) < epsilon);
396 
397  if (nbits > 0)
398  nbits_ = nbits;
399  else
400  nbits_ = p1->nbits() - lr_;
401 
402  Kmap_["2"] = Kmap_["2"] + lr_;
403  }
404 
405  ~VarAdjustKR() override = default;
406 
407  void local_calculate() override;
408  void print(std::ofstream &fs, Verilog, int l1 = 0, int l2 = 0, int l3 = 0) override;
409  void print(std::ofstream &fs, HLS, int l1 = 0, int l2 = 0, int l3 = 0) override;
410 
411  protected:
412  int lr_;
413  };
414 
415  class VarParam : public VarBase {
416  public:
417  VarParam(imathGlobals *globals, std::string name, double fval, int nbits)
418  : VarBase(globals, name, nullptr, nullptr, nullptr, 0) {
419  op_ = "const";
420  nbits_ = nbits;
421  int l = log2(std::abs(fval)) + 1.9999999 - nbits;
422  Kmap_["2"] = l;
423  K_ = pow(2, l);
424  fval_ = fval;
425  ival_ = fval / K_;
426  }
428  : VarBase(globals, name, nullptr, nullptr, nullptr, 0) {
429  op_ = "const";
430  K_ = K;
431  nbits_ = log2(fval / K) + 1.999999; //plus one to round up
432  if (!units.empty())
433  Kmap_[units] = 1;
434  else {
435  //defining a constant, K should be a power of two
436  int l = log2(K);
437  if (std::abs(pow(2, l) / K - 1) > 1e-5) {
438  char slog[100];
439  snprintf(slog, 100, "defining unitless constant, yet K is not a power of 2! %g, %g", K, pow(2, l));
440  edm::LogVerbatim("Tracklet") << slog;
441  }
442  Kmap_["2"] = l;
443  }
444  }
445 
446  ~VarParam() override = default;
447 
448  void set_fval(double fval) {
449  fval_ = fval;
450  if (fval > 0)
451  ival_ = fval / K_ + 0.5;
452  else
453  ival_ = fval / K_ - 0.5;
454  val_ = ival_ * K_;
455  }
456  void set_ival(int ival) {
457  ival_ = ival;
458  fval_ = ival * K_;
459  val_ = fval_;
460  }
461  void print(std::ofstream &fs, Verilog, int l1 = 0, int l2 = 0, int l3 = 0) override;
462  void print(std::ofstream &fs, HLS, int l1 = 0, int l2 = 0, int l3 = 0) override;
463  };
464 
465  class VarDef : public VarBase {
466  public:
467  //construct from scratch
468  VarDef(imathGlobals *globals, std::string name, std::string units, double fmax, double K)
469  : VarBase(globals, name, nullptr, nullptr, nullptr, 1) {
470  op_ = "def";
471  K_ = K;
472  nbits_ = log2(fmax / K) + 1.999999; //plus one to round up
473  if (!units.empty())
474  Kmap_[units] = 1;
475  else {
476  //defining a constant, K should be a power of two
477  int l = log2(K);
478  if (std::abs(pow(2, l) / K - 1) > 1e-5) {
479  char slog[100];
480  snprintf(slog, 100, "defining unitless constant, yet K is not a power of 2! %g, %g", K, pow(2, l));
481  edm::LogVerbatim("Tracklet") << slog;
482  }
483  Kmap_["2"] = l;
484  }
485  }
486  //construct from abother variable (all provenance info is lost!)
487  VarDef(imathGlobals *globals, std::string name, VarBase *p) : VarBase(globals, name, nullptr, nullptr, nullptr, 1) {
488  op_ = "def";
489  K_ = p->K();
490  nbits_ = p->nbits();
491  Kmap_ = p->Kmap();
492  }
493  void set_fval(double fval) {
494  fval_ = fval;
495  if (fval > 0)
496  ival_ = fval / K_;
497  else
498  ival_ = fval / K_ - 1;
499  val_ = ival_ * K_;
500  }
501  void set_ival(int ival) {
502  ival_ = ival;
503  fval_ = ival * K_;
504  val_ = ival_ * K_;
505  }
506  ~VarDef() override = default;
507  void print(std::ofstream &fs, Verilog, int l1 = 0, int l2 = 0, int l3 = 0) override;
508  void print(std::ofstream &fs, HLS, int l1 = 0, int l2 = 0, int l3 = 0) override;
509  };
510 
511  class VarAdd : public VarBase {
512  public:
513  VarAdd(imathGlobals *globals, std::string name, VarBase *p1, VarBase *p2, double range = -1, int nmax = 18)
514  : VarBase(globals, name, p1, p2, nullptr, 1) {
515  op_ = "add";
516 
517  std::map<std::string, int> map1 = p1->Kmap();
518  std::map<std::string, int> map2 = p2->Kmap();
519  int s1 = map1["2"];
520  int s2 = map2["2"];
521 
522  //first check if the constants are all lined up
523  //go over the two maps subtracting the units
524  for (const auto &it : map2) {
525  if (map1.find(it.first) == map1.end())
526  map1[it.first] = -it.second;
527  else
528  map1[it.first] = map1[it.first] - it.second;
529  }
530 
531  char slog[100];
532 
533  //assert if different
534  for (const auto &it : map1) {
535  if (it.second != 0) {
536  if (it.first != "2") {
537  snprintf(
538  slog, 100, "VarAdd: bad units! %s^%i for variable %s", (it.first).c_str(), it.second, name_.c_str());
539  edm::LogVerbatim("Tracklet") << slog;
540  p1->dump_msg();
541  p2->dump_msg();
542  throw cms::Exception("BadConfig") << "imath units are different!";
543  }
544  }
545  }
546 
547  double ki1 = p1->K() / pow(2, s1);
548  double ki2 = p2->K() / pow(2, s2);
549  //those should be the same
550  if (std::abs(ki1 / ki2 - 1.) > 1e-6) {
551  snprintf(slog, 100, "VarAdd: bad constants! %f %f for variable %s", ki1, ki2, name_.c_str());
552  edm::LogVerbatim("Tracklet") << slog;
553  p1->dump_msg();
554  p2->dump_msg();
555  throw cms::Exception("BadConfig") << "imath constants are different!";
556  }
557  //everything checks out!
558 
559  Kmap_ = p1->Kmap();
560 
561  int s0 = s1 < s2 ? s1 : s2;
562  shift1 = s1 - s0;
563  shift2 = s2 - s0;
564 
565  int n1 = p1->nbits() + shift1;
566  int n2 = p2->nbits() + shift2;
567  int n0 = 1 + (n1 > n2 ? n1 : n2);
568 
569  //before shifting, check the range
570  if (range > 0) {
571  n0 = log2(range / ki1 / pow(2, s0)) + 1e-9;
572  n0 = n0 + 2;
573  }
574 
575  if (n0 <= nmax) { //if it fits, we're done
576  ps_ = 0;
577  Kmap_["2"] = s0;
578  nbits_ = n0;
579  } else {
580  ps_ = n0 - nmax;
581  Kmap_["2"] = s0 + ps_;
582  nbits_ = nmax;
583  }
584 
585  K_ = ki1 * pow(2, Kmap_["2"]);
586  }
587  ~VarAdd() override = default;
588  void local_calculate() override;
589  void print(std::ofstream &fs, Verilog, int l1 = 0, int l2 = 0, int l3 = 0) override;
590  void print(std::ofstream &fs, HLS, int l1 = 0, int l2 = 0, int l3 = 0) override;
591 
592  protected:
593  int ps_;
594  int shift1;
595  int shift2;
596  };
597 
598  class VarSubtract : public VarBase {
599  public:
600  VarSubtract(imathGlobals *globals, std::string name, VarBase *p1, VarBase *p2, double range = -1, int nmax = 18)
601  : VarBase(globals, name, p1, p2, nullptr, 1) {
602  op_ = "subtract";
603 
604  std::map<std::string, int> map1 = p1->Kmap();
605  std::map<std::string, int> map2 = p2->Kmap();
606  int s1 = map1["2"];
607  int s2 = map2["2"];
608 
609  //first check if the constants are all lined up go over the two maps subtracting the units
610  for (const auto &it : map2) {
611  if (map1.find(it.first) == map1.end())
612  map1[it.first] = -it.second;
613  else
614  map1[it.first] = map1[it.first] - it.second;
615  }
616 
617  char slog[100];
618 
619  //assert if different
620  for (const auto &it : map1) {
621  if (it.second != 0) {
622  if (it.first != "2") {
623  snprintf(
624  slog, 100, "VarAdd: bad units! %s^%i for variable %s", (it.first).c_str(), it.second, name_.c_str());
625  edm::LogVerbatim("Tracklet") << slog;
626  p1->dump_msg();
627  p2->dump_msg();
628  throw cms::Exception("BadConfig") << "imath units are different!";
629  }
630  }
631  }
632 
633  double ki1 = p1->K() / pow(2, s1);
634  double ki2 = p2->K() / pow(2, s2);
635  //those should be the same
636  if (std::abs(ki1 / ki2 - 1.) > 1e-6) {
637  snprintf(slog, 100, "VarAdd: bad constants! %f %f for variable %s", ki1, ki2, name_.c_str());
638  edm::LogVerbatim("Tracklet") << slog;
639  p1->dump_msg();
640  p2->dump_msg();
641  throw cms::Exception("BadConfig") << "imath constants are different!";
642  }
643  //everything checks out!
644 
645  Kmap_ = p1->Kmap();
646 
647  int s0 = s1 < s2 ? s1 : s2;
648  shift1 = s1 - s0;
649  shift2 = s2 - s0;
650 
651  int n1 = p1->nbits() + shift1;
652  int n2 = p2->nbits() + shift2;
653  int n0 = 1 + (n1 > n2 ? n1 : n2);
654 
655  //before shifting, check the range
656  if (range > 0) {
657  n0 = log2(range / ki1 / pow(2, s0)) + 1e-9;
658  n0 = n0 + 2;
659  }
660 
661  if (n0 <= nmax) { //if it fits, we're done
662  ps_ = 0;
663  Kmap_["2"] = s0;
664  nbits_ = n0;
665  } else {
666  ps_ = n0 - nmax;
667  Kmap_["2"] = s0 + ps_;
668  nbits_ = nmax;
669  }
670 
671  K_ = ki1 * pow(2, Kmap_["2"]);
672  }
673 
674  ~VarSubtract() override = default;
675 
676  void local_calculate() override;
677  void print(std::ofstream &fs, Verilog, int l1 = 0, int l2 = 0, int l3 = 0) override;
678  void print(std::ofstream &fs, HLS, int l1 = 0, int l2 = 0, int l3 = 0) override;
679 
680  protected:
681  int ps_;
682  int shift1;
683  int shift2;
684  };
685 
686  class VarNounits : public VarBase {
687  public:
688  VarNounits(imathGlobals *globals, std::string name, VarBase *p1, int ps = 17)
689  : VarBase(globals, name, p1, nullptr, nullptr, MULT_LATENCY) {
690  op_ = "nounits";
691  ps_ = ps;
692  nbits_ = p1->nbits();
693 
694  int s1 = p1->shift();
695  double ki = p1->K() / pow(2, s1);
696  int m = log2(ki);
697 
698  K_ = pow(2, s1 + m);
699  Kmap_["2"] = s1 + m;
700  double c = ki * pow(2, -m);
701  cI_ = c * pow(2, ps_);
702  }
703  ~VarNounits() override = default;
704 
705  void local_calculate() override;
706  void print(std::ofstream &fs, Verilog, int l1 = 0, int l2 = 0, int l3 = 0) override;
707  void print(std::ofstream &fs, HLS, int l1 = 0, int l2 = 0, int l3 = 0) override;
708 
709  protected:
710  int ps_;
711  int cI_;
712  };
713 
714  class VarShiftround : public VarBase {
715  public:
717  : VarBase(globals, name, p1, nullptr, nullptr, 1) { // latency is one because there is an addition
718  op_ = "shiftround";
719  shift_ = shift;
720 
721  nbits_ = p1->nbits() - shift;
722  Kmap_ = p1->Kmap();
723  K_ = p1->K();
724  }
725  ~VarShiftround() override = default;
726 
727  void local_calculate() override;
728  void print(std::ofstream &fs, Verilog, int l1 = 0, int l2 = 0, int l3 = 0) override;
729  void print(std::ofstream &fs, HLS, int l1 = 0, int l2 = 0, int l3 = 0) override;
730 
731  protected:
732  int shift_;
733  };
734 
735  class VarShift : public VarBase {
736  public:
738  : VarBase(globals, name, p1, nullptr, nullptr, 0) {
739  op_ = "shift";
740  shift_ = shift;
741 
742  nbits_ = p1->nbits() - shift;
743  Kmap_ = p1->Kmap();
744  K_ = p1->K();
745  }
746  ~VarShift() override = default;
747  void local_calculate() override;
748  void print(std::ofstream &fs, Verilog, int l1 = 0, int l2 = 0, int l3 = 0) override;
749  void print(std::ofstream &fs, HLS, int l1 = 0, int l2 = 0, int l3 = 0) override;
750 
751  protected:
752  int shift_;
753  };
754 
755  class VarNeg : public VarBase {
756  public:
757  VarNeg(imathGlobals *globals, std::string name, VarBase *p1) : VarBase(globals, name, p1, nullptr, nullptr, 1) {
758  op_ = "neg";
759  nbits_ = p1->nbits();
760  Kmap_ = p1->Kmap();
761  K_ = p1->K();
762  }
763  ~VarNeg() override = default;
764  void local_calculate() override;
765  void print(std::ofstream &fs, Verilog, int l1 = 0, int l2 = 0, int l3 = 0) override;
766  void print(std::ofstream &fs, HLS, int l1 = 0, int l2 = 0, int l3 = 0) override;
767  };
768 
769  class VarTimesC : public VarBase {
770  public:
771  VarTimesC(imathGlobals *globals, std::string name, VarBase *p1, double cF, int ps = 17)
772  : VarBase(globals, name, p1, nullptr, nullptr, MULT_LATENCY) {
773  op_ = "timesC";
774  cF_ = cF;
775  ps_ = ps;
776 
777  nbits_ = p1->nbits();
778  Kmap_ = p1->Kmap();
779  K_ = p1->K();
780 
781  int s1 = Kmap_["2"];
782  double l = log2(std::abs(cF));
783  if (l > 0)
784  l += 0.999999;
785  int m = l;
786 
787  cI_ = cF * pow(2, ps - m);
788  K_ = K_ * pow(2, m);
789  Kmap_["2"] = s1 + m;
790  }
791  ~VarTimesC() override = default;
792  void local_calculate() override;
793  void print(std::ofstream &fs, Verilog, int l1 = 0, int l2 = 0, int l3 = 0) override;
794  void print(std::ofstream &fs, HLS, int l1 = 0, int l2 = 0, int l3 = 0) override;
795 
796  protected:
797  int ps_;
798  int cI_;
799  double cF_;
800  };
801 
802  class VarMult : public VarBase {
803  public:
804  VarMult(imathGlobals *globals, std::string name, VarBase *p1, VarBase *p2, double range = -1, int nmax = 18)
805  : VarBase(globals, name, p1, p2, nullptr, MULT_LATENCY) {
806  op_ = "mult";
807 
808  const std::map<std::string, int> map1 = p1->Kmap();
809  const std::map<std::string, int> map2 = p2->Kmap();
810  for (const auto &it : map1) {
811  if (Kmap_.find(it.first) == Kmap_.end())
812  Kmap_[it.first] = it.second;
813  else
814  Kmap_[it.first] = Kmap_[it.first] + it.second;
815  }
816  for (const auto &it : map2) {
817  if (Kmap_.find(it.first) == Kmap_.end())
818  Kmap_[it.first] = it.second;
819  else
820  Kmap_[it.first] = Kmap_[it.first] + it.second;
821  }
822  K_ = p1->K() * p2->K();
823  int s0 = Kmap_["2"];
824 
825  int n0 = p1->nbits() + p2->nbits();
826  if (range > 0) {
827  n0 = log2(range / K_) + 1e-9;
828  n0 = n0 + 2;
829  }
830  if (n0 < nmax) {
831  ps_ = 0;
832  nbits_ = n0;
833  } else {
834  ps_ = n0 - nmax;
835  nbits_ = nmax;
836  Kmap_["2"] = s0 + ps_;
837  K_ = K_ * pow(2, ps_);
838  }
839  }
840  ~VarMult() override = default;
841  void local_calculate() override;
842  void print(std::ofstream &fs, Verilog, int l1 = 0, int l2 = 0, int l3 = 0) override;
843  void print(std::ofstream &fs, HLS, int l1 = 0, int l2 = 0, int l3 = 0) override;
844 
845  protected:
846  int ps_;
847  };
848 
849  class VarDSPPostadd : public VarBase {
850  public:
852  imathGlobals *globals, std::string name, VarBase *p1, VarBase *p2, VarBase *p3, double range = -1, int nmax = 18)
853  : VarBase(globals, name, p1, p2, p3, DSP_LATENCY) {
854  op_ = "DSP_postadd";
855 
856  //first, get constants for the p1*p2
857  std::map<std::string, int> map1 = p1->Kmap();
858  std::map<std::string, int> map2 = p2->Kmap();
859  for (const auto &it : map2) {
860  if (map1.find(it.first) == map1.end())
861  map1[it.first] = it.second;
862  else
863  map1[it.first] = map1[it.first] + it.second;
864  }
865  double k0 = p1->K() * p2->K();
866  int s0 = map1["2"];
867 
868  //now addition
869  std::map<std::string, int> map3 = p3->Kmap();
870  int s3 = map3["2"];
871 
872  //first check if the constants are all lined up
873  //go over the two maps subtracting the units
874  for (const auto &it : map3) {
875  if (map1.find(it.first) == map1.end())
876  map1[it.first] = -it.second;
877  else
878  map1[it.first] = map1[it.first] - it.second;
879  }
880 
881  char slog[100];
882 
883  //assert if different
884  for (const auto &it : map1) {
885  if (it.second != 0) {
886  if (it.first != "2") {
887  snprintf(slog,
888  100,
889  "VarDSPPostadd: bad units! %s^%i for variable %s",
890  (it.first).c_str(),
891  it.second,
892  name_.c_str());
893  edm::LogVerbatim("Tracklet") << slog;
894  p1->dump_msg();
895  p2->dump_msg();
896  p3->dump_msg();
897  throw cms::Exception("BadConfig") << "imath units are different!";
898  }
899  }
900  }
901 
902  double ki1 = k0 / pow(2, s0);
903  double ki2 = p3->K() / pow(2, s3);
904  //those should be the same
905  if (std::abs(ki1 / ki2 - 1.) > 1e-6) {
906  snprintf(slog, 100, "VarDSPPostadd: bad constants! %f %f for variable %s", ki1, ki2, name_.c_str());
907  edm::LogVerbatim("Tracklet") << slog;
908  p1->dump_msg();
909  p2->dump_msg();
910  p3->dump_msg();
911  throw cms::Exception("BadConfig") << "imath constants are different!";
912  }
913  //everything checks out!
914 
915  shift3_ = s3 - s0;
916  if (shift3_ < 0) {
917  throw cms::Exception("BadConfig") << "imath VarDSPPostadd: loosing precision on C in A*B+C: " << shift3_;
918  }
919 
920  Kmap_ = p3->Kmap();
921  Kmap_["2"] = Kmap_["2"] - shift3_;
922 
923  int n12 = p1->nbits() + p2->nbits();
924  int n3 = p3->nbits() + shift3_;
925  int n0 = 1 + (n12 > n3 ? n12 : n3);
926 
927  //before shifting, check the range
928  if (range > 0) {
929  n0 = log2(range / ki2 / pow(2, s3)) + 1e-9;
930  n0 = n0 + 2;
931  }
932 
933  if (n0 <= nmax) { //if it fits, we're done
934  ps_ = 0;
935  nbits_ = n0;
936  } else {
937  ps_ = n0 - nmax;
938  Kmap_["2"] = Kmap_["2"] + ps_;
939  nbits_ = nmax;
940  }
941 
942  K_ = ki2 * pow(2, Kmap_["2"]);
943  }
944  ~VarDSPPostadd() override = default;
945 
946  void local_calculate() override;
947  using VarBase::print;
948  void print(std::ofstream &fs, Verilog, int l1 = 0, int l2 = 0, int l3 = 0) override;
949 
950  protected:
951  int ps_;
952  int shift3_;
953  };
954 
955  class VarInv : public VarBase {
956  public:
957  enum mode { pos, neg, both };
958 
961  VarBase *p1,
962  double offset,
963  int nbits,
964  int n,
965  unsigned int shift,
966  mode m,
967  int nbaddr = -1)
968  : VarBase(globals, name, p1, nullptr, nullptr, LUT_LATENCY) {
969  op_ = "inv";
970  offset_ = offset;
971  nbits_ = nbits;
972  n_ = n;
973  shift_ = shift;
974  m_ = m;
975  if (nbaddr < 0)
976  nbaddr = p1->nbits();
977  nbaddr_ = nbaddr - shift;
978  if (m_ != mode::both)
979  nbaddr_--;
980  Nelements_ = 1 << nbaddr_;
981  mask_ = Nelements_ - 1;
982  ashift_ = sizeof(int) * 8 - nbaddr_;
983 
984  const std::map<std::string, int> map1 = p1->Kmap();
985  for (const auto &it : map1)
986  Kmap_[it.first] = -it.second;
987  Kmap_["2"] = Kmap_["2"] - n;
988  K_ = pow(2, -n) / p1->K();
989 
990  LUT = new int[Nelements_];
991  double offsetI = lround(offset_ / p1_->K());
992  for (int i = 0; i < Nelements_; ++i) {
993  int i1 = addr_to_ival(i);
994  LUT[i] = gen_inv(offsetI + i1);
995  }
996  }
997  ~VarInv() override { delete[] LUT; }
998 
999  void set_mode(mode m) { m_ = m; }
1000  void initLUT(double offset);
1001  double offset() { return offset_; }
1002  double Ioffset() { return offset_ / p1_->K(); }
1003 
1004  void local_calculate() override;
1005  void print(std::ofstream &fs, Verilog, int l1 = 0, int l2 = 0, int l3 = 0) override;
1006  void print(std::ofstream &fs, HLS, int l1 = 0, int l2 = 0, int l3 = 0) override;
1007  void writeLUT(std::ofstream &fs) const { writeLUT(fs, verilog); }
1008  void writeLUT(std::ofstream &fs, Verilog) const;
1009  void writeLUT(std::ofstream &fs, HLS) const;
1010 
1011  int ival_to_addr(int ival) { return ((ival >> shift_) & mask_); }
1012  int addr_to_ival(int addr) {
1013  switch (m_) {
1014  case mode::pos:
1015  return addr << shift_;
1016  case mode::neg:
1017  return (addr - Nelements_) << shift_;
1018  case mode::both:
1019  return (addr << ashift_) >> (ashift_ - shift_);
1020  }
1021  assert(0);
1022  }
1023  int gen_inv(int i) {
1024  unsigned int ms = sizeof(int) * 8 - nbits_;
1025  int lut = 0;
1026  if (i > 0) {
1027  int i1 = i + (1 << shift_) - 1;
1028  int lut1 = (lround((1 << n_) / i) << ms) >> ms;
1029  int lut2 = (lround((1 << n_) / (i1)) << ms) >> ms;
1030  lut = 0.5 * (lut1 + lut2);
1031  } else if (i < -1) {
1032  int i1 = i + (1 << shift_) - 1;
1033  int i2 = i;
1034  int lut1 = (lround((1 << n_) / i1) << ms) >> ms;
1035  int lut2 = (lround((1 << n_) / i2) << ms) >> ms;
1036  lut = 0.5 * (lut1 + lut2);
1037  }
1038  return lut;
1039  }
1040 
1041  protected:
1042  double offset_;
1043  int n_;
1045  unsigned int shift_;
1046  unsigned int mask_;
1047  unsigned int ashift_;
1049  int nbaddr_;
1050 
1051  int *LUT;
1052  };
1053 
1054  class VarCut : public VarBase {
1055  public:
1056  VarCut(imathGlobals *globals, double lower_cut, double upper_cut)
1057  : VarBase(globals, "", nullptr, nullptr, nullptr, 0),
1060  parent_flag_(nullptr) {
1061  op_ = "cut";
1062  }
1063 
1065  : VarCut(globals, lower_cut, upper_cut) {
1067  }
1068  ~VarCut() override = default;
1069 
1070  double lower_cut() const { return lower_cut_; }
1071  double upper_cut() const { return upper_cut_; }
1072 
1073  void local_passes(std::map<const VarBase *, std::vector<bool> > &passes,
1074  const std::map<const VarBase *, std::vector<bool> > *const previous_passes = nullptr) const;
1075  using VarBase::print;
1076  void print(std::map<const VarBase *, std::set<std::string> > &cut_strings,
1077  const int step,
1078  Verilog,
1079  const std::map<const VarBase *, std::set<std::string> > *const previous_cut_strings = nullptr) const;
1080  void print(std::map<const VarBase *, std::set<std::string> > &cut_strings,
1081  const int step,
1082  HLS,
1083  const std::map<const VarBase *, std::set<std::string> > *const previous_cut_strings = nullptr) const;
1084 
1085  void set_parent_flag(VarFlag *parent_flag, const bool call_add_cut);
1087  void set_cut_var(VarBase *cut_var, const bool call_add_cut = true);
1088 
1089  protected:
1090  double lower_cut_;
1091  double upper_cut_;
1093  };
1094 
1095  class VarFlag : public VarBase {
1096  public:
1097  template <class... Args>
1099  : VarBase(globals, name, nullptr, nullptr, nullptr, 0) {
1100  op_ = "flag";
1101  nbits_ = 1;
1102  add_cuts(cut, args...);
1103  }
1104 
1105  template <class... Args>
1106  void add_cuts(VarBase *cut, Args... args) {
1107  add_cut(cut);
1108  add_cuts(args...);
1109  }
1110 
1112 
1113  void add_cut(VarBase *cut, const bool call_set_parent_flag = true);
1114 
1115  void calculate_step();
1116  bool passes();
1117  void print(std::ofstream &fs, Verilog, int l1 = 0, int l2 = 0, int l3 = 0) override;
1118  void print(std::ofstream &fs, HLS, int l1 = 0, int l2 = 0, int l3 = 0) override;
1119  };
1120 }; // namespace trklet
1121 #endif
void set_cut_var(VarBase *cut_var, const bool call_add_cut=true)
Definition: imath.cc:365
VarShiftround(imathGlobals *globals, std::string name, VarBase *p1, int shift)
Definition: imath.h:716
int * LUT
Definition: imath.h:1051
Log< level::Info, true > LogVerbatim
std::string name() const
Definition: imath.h:207
int latency_
Definition: imath.h:301
void print(std::ofstream &fs, Verilog, int l1=0, int l2=0, int l3=0) override
void print(std::ofstream &fs, Verilog, int l1=0, int l2=0, int l3=0) override
static void verilog_print(const std::vector< VarBase *> &v, std::ofstream &fs)
Definition: imath.h:274
double fval() const
Definition: imath.h:212
static void hls_print(const std::vector< VarBase *> &v, std::ofstream &fs)
Definition: imath.h:275
int gen_inv(int i)
Definition: imath.h:1023
void print(std::ofstream &fs, Verilog, int l1=0, int l2=0, int l3=0) override
void pipe_increment()
Definition: imath.h:271
void analyze()
Definition: imath.cc:29
VarFlag * parent_flag_
Definition: imath.h:1092
VarDef(imathGlobals *globals, std::string name, VarBase *p)
Definition: imath.h:487
~VarMult() override=default
void inputs(std::vector< VarBase *> *vd)
Definition: imath.cc:205
int pipe_counter()
Definition: imath.h:270
void add_cut(VarCut *cut, const bool call_set_cut_var=true)
Definition: imath.cc:359
unsigned int mask_
Definition: imath.h:1046
void local_calculate() override
double val_
Definition: imath.h:306
void print_step(int step, std::ofstream &fs, Verilog)
void local_calculate() override
void add_cuts(VarBase *cut, Args... args)
Definition: imath.h:1106
void print_all(std::ofstream &fs, Verilog)
double range() const
Definition: imath.h:245
int nbits() const
Definition: imath.h:243
VarBase * cut_var()
Definition: imath.cc:388
void print(std::ofstream &fs, Verilog, int l1=0, int l2=0, int l3=0) override
VarFlag(imathGlobals *globals, std::string name, VarBase *cut, Args... args)
Definition: imath.h:1098
~VarTimesC() override=default
VarBase * p1_
Definition: imath.h:297
double offset()
Definition: imath.h:1001
void local_calculate() override
void add_cuts(VarBase *cut)
Definition: imath.h:1111
std::string pipe_delays(const int step)
Definition: imath.cc:184
bool passes()
Definition: imath.cc:395
VarBase * cut_var_
Definition: imath.h:309
double Ioffset()
Definition: imath.h:1002
imathGlobals * globals_
Definition: imath.h:295
VarBase * p2_
Definition: imath.h:298
bool has_delay(int i)
Definition: imath.cc:165
void set_fval(double fval)
Definition: imath.h:493
std::string kstring() const
Definition: imath.cc:18
VarDef(imathGlobals *globals, std::string name, std::string units, double fmax, double K)
Definition: imath.h:468
void print(std::ofstream &fs, Verilog, int l1=0, int l2=0, int l3=0) override
assert(be >=bs)
~VarAdd() override=default
~VarParam() override=default
std::vector< int > pipe_delays_
Definition: imath.h:316
void print(std::map< const VarBase *, std::set< std::string > > &cut_strings, const int step, Verilog, const std::map< const VarBase *, std::set< std::string > > *const previous_cut_strings=nullptr) const
static struct trklet::VarBase::Verilog verilog
VarInv(imathGlobals *globals, std::string name, VarBase *p1, double offset, int nbits, int n, unsigned int shift, mode m, int nbaddr=-1)
Definition: imath.h:959
int shift() const
Definition: imath.h:247
int ival_to_addr(int ival)
Definition: imath.h:1011
void passes(std::map< const VarBase *, std::vector< bool > > &passes, const std::map< const VarBase *, std::vector< bool > > *const previous_passes=nullptr) const
Definition: imath.cc:332
int addr_to_ival(int addr)
Definition: imath.h:1012
VarBase(imathGlobals *globals, std::string name, VarBase *p1, VarBase *p2, VarBase *p3, int l)
Definition: imath.h:161
void local_calculate() override
VarBase * p3_
Definition: imath.h:299
VarSubtract(imathGlobals *globals, std::string name, VarBase *p1, VarBase *p2, double range=-1, int nmax=18)
Definition: imath.h:600
std::string op_
Definition: imath.h:300
int Nelements_
Definition: imath.h:1048
VarShift(imathGlobals *globals, std::string name, VarBase *p1, int shift)
Definition: imath.h:737
VarAdd(imathGlobals *globals, std::string name, VarBase *p1, VarBase *p2, double range=-1, int nmax=18)
Definition: imath.h:513
void print(std::ofstream &fs, Verilog, int l1=0, int l2=0, int l3=0) override
void initLUT(double offset)
Definition: imath.cc:142
VarParam(imathGlobals *globals, std::string name, std::string units, double fval, double K)
Definition: imath.h:427
double offset_
Definition: imath.h:1042
void local_calculate() override
#define MULT_LATENCY
Definition: imath.h:139
virtual void local_calculate()
Definition: imath.h:255
double lower_cut_
Definition: imath.h:1090
bool local_passes() const
Definition: imath.cc:315
~VarShift() override=default
int n0
Definition: AMPTWrapper.h:44
void set_fval(double fval)
Definition: imath.h:448
unsigned int ashift_
Definition: imath.h:1047
VarBase * p1() const
Definition: imath.h:209
std::map< std::string, int > Kmap() const
Definition: imath.h:244
void local_calculate() override
#define LUT_LATENCY
Definition: imath.h:140
void local_calculate() override
~VarSubtract() override=default
int step() const
Definition: imath.h:250
void print(std::ofstream &fs, Verilog, int l1=0, int l2=0, int l3=0) override
std::string op() const
Definition: imath.h:208
~VarShiftround() override=default
int latency() const
Definition: imath.h:251
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void print(std::ofstream &fs, Verilog, int l1=0, int l2=0, int l3=0) override
void dump_msg()
Definition: imath.cc:96
VarAdjustK(imathGlobals *globals, std::string name, VarBase *p1, double Knew, double epsilon=1e-5, bool do_assert=false, int nbits=-1)
Definition: imath.h:337
VarNounits(imathGlobals *globals, std::string name, VarBase *p1, int ps=17)
Definition: imath.h:688
void local_calculate() override
double maxval() const
Definition: imath.h:230
static std::string itos(int i)
Definition: imath.cc:16
void adjust(double Knew, double epsilon=1e-5, bool do_assert=false, int nbits=-1)
Definition: imath.cc:120
long int ival_
Definition: imath.h:305
VarBase * p3() const
Definition: imath.h:211
std::vector< DeviationSensor2D * > vd
double upper_cut() const
Definition: imath.h:1071
void set_parent_flag(VarFlag *parent_flag, const bool call_add_cut)
Definition: imath.cc:382
#define DSP_LATENCY
Definition: imath.h:141
void local_calculate() override
virtual void print(std::ofstream &fs, HLS, int l1=0, int l2=0, int l3=0)
Definition: imath.h:259
double K() const
Definition: imath.h:246
void local_calculate() override
double minval_
Definition: imath.h:322
void print(std::ofstream &fs, Verilog, int l1=0, int l2=0, int l3=0) override
void print(std::ofstream &fs, Verilog, int l1=0, int l2=0, int l3=0) override
virtual void print(std::ofstream &fs, Verilog, int l1=0, int l2=0, int l3=0)
Definition: imath.h:256
static struct trklet::VarBase::HLS hls
bool readytoprint_
Definition: imath.h:319
void add_latency(unsigned int l)
Definition: imath.h:252
void set_mode(mode m)
Definition: imath.h:999
void print(std::ofstream &fs, Verilog, int l1=0, int l2=0, int l3=0) override
~VarDSPPostadd() override=default
void add_cut(VarBase *cut, const bool call_set_parent_flag=true)
Definition: imath.cc:373
VarDSPPostadd(imathGlobals *globals, std::string name, VarBase *p1, VarBase *p2, VarBase *p3, double range=-1, int nmax=18)
Definition: imath.h:851
void print(std::ofstream &fs, Verilog, int l1=0, int l2=0, int l3=0) override
void print_truncation(std::string &t, const std::string &o1, const int ps, Verilog) const
int pipe_counter_
Definition: imath.h:315
VarAdjustKR(imathGlobals *globals, std::string name, VarBase *p1, double Knew, double epsilon=1e-5, bool do_assert=false, int nbits=-1)
Definition: imath.h:378
double upper_cut_
Definition: imath.h:1091
void makeready()
Definition: imath.cc:151
static std::string pipe_delay_wire(VarBase *v, std::string name_delayed, int nbits, int delay)
Definition: imath.cc:196
long int ival() const
Definition: imath.h:213
std::vector< VarBase * > cuts_
Definition: imath.h:308
std::string name_
Definition: imath.h:296
static std::string pipe_delay(VarBase *v, int nbits, int delay)
Definition: imath.cc:173
std::map< std::string, int > Kmap_
Definition: imath.h:313
VarBase * p2() const
Definition: imath.h:210
static void design_print(const std::vector< VarBase *> &v, std::ofstream &fs, Verilog)
VarCut(imathGlobals *globals, double lower_cut, double upper_cut)
Definition: imath.h:1056
void reset()
Definition: imath.h:235
~VarCut() override=default
double maxval_
Definition: imath.h:323
void local_calculate() override
~VarAdjustK() override=default
~VarNeg() override=default
VarTimesC(imathGlobals *globals, std::string name, VarBase *p1, double cF, int ps=17)
Definition: imath.h:771
bool calculate()
Definition: imath.h:254
double lower_cut() const
Definition: imath.h:1070
void set_ival(int ival)
Definition: imath.h:456
TString units(TString variable, Char_t axis)
VarCut(imathGlobals *globals, VarBase *cut_var, double lower_cut, double upper_cut)
Definition: imath.h:1064
double K_
Definition: imath.h:312
double cF_
Definition: imath.h:799
void set_ival(int ival)
Definition: imath.h:501
void print(std::ofstream &fs, Verilog, int l1=0, int l2=0, int l3=0) override
int shift1
Definition: imath.h:594
double minval() const
Definition: imath.h:229
VarMult(imathGlobals *globals, std::string name, VarBase *p1, VarBase *p2, double range=-1, int nmax=18)
Definition: imath.h:804
~VarAdjustKR() override=default
~VarInv() override
Definition: imath.h:997
int shift2
Definition: imath.h:595
bool readytoanalyze_
Definition: imath.h:318
VarNeg(imathGlobals *globals, std::string name, VarBase *p1)
Definition: imath.h:757
~VarDef() override=default
virtual ~VarBase()
Definition: imath.h:191
VarParam(imathGlobals *globals, std::string name, double fval, int nbits)
Definition: imath.h:417
void print(std::ofstream &fs, Verilog, int l1=0, int l2=0, int l3=0) override
std::string dump()
Definition: imath.cc:78
step
Definition: StallMonitor.cc:98
bool usedasinput_
Definition: imath.h:320
void writeLUT(std::ofstream &fs) const
Definition: imath.h:1007
std::vector< unsigned short int > LUT
Definition: DTTracoLUTs.h:31
unsigned int shift_
Definition: imath.h:1045
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
double fval_
Definition: imath.h:304
VarFlag * parent_flag()
Definition: imath.h:1086
void add_delay(int i)
Definition: imath.h:272
~VarNounits() override=default
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
void print(std::ofstream &fs, Verilog, int l1=0, int l2=0, int l3=0) override
void print_cuts(std::map< const VarBase *, std::set< std::string > > &cut_strings, const int step, Verilog, const std::map< const VarBase *, std::set< std::string > > *const previous_cut_strings=nullptr) const
void local_calculate() override