CORSIKA  @c8_version@
The framework to simulate particle cascades for astroparticle physics
quantity.hpp
Go to the documentation of this file.
1 
21 /*
22  * Unless otherwise specified, the definitions of all units in this
23  * file are from NIST Special Publication 811, found online at
24  * http://physics.nist.gov/Document/sp811.pdf
25  * Other sources: OED = Oxford English Dictionary
26  */
27 
28 #ifndef PHYS_UNITS_QUANTITY_HPP_INCLUDED
29 #define PHYS_UNITS_QUANTITY_HPP_INCLUDED
30 
31 #include <cmath>
32 #include <cstdlib>
33 #include <utility> // std::declval
34 #include <type_traits>
35 
37 
38 namespace phys {
39 
41 
42  namespace units {
43 
44 #ifdef PHYS_UNITS_REP_TYPE
45  using Rep = PHYS_UNITS_REP_TYPE;
46 #else
47  using Rep = double;
48 #endif
49 
50  /*
51  * declare now, define later.
52  */
53  template <typename Dims, typename T = Rep>
54  class quantity;
55 
60  template <int D1, int D2, int D3, int D4 = 0, int D5 = 0, int D6 = 0, int D7 = 0,
61  int D8 = 0>
62  struct dimensions {
63  enum {
64  dim1 = D1,
65  dim2 = D2,
66  dim3 = D3,
67  dim4 = D4,
68  dim5 = D5,
69  dim6 = D6,
70  dim7 = D7,
71  dim8 = D8,
72 
73  is_all_zero = D1 == 0 && D2 == 0 && D3 == 0 && D4 == 0 && D5 == 0 && D6 == 0 &&
74  D7 == 0 && D8 == 0,
75 
76  is_base = 1 == (D1 != 0) + (D2 != 0) + (D3 != 0) + (D4 != 0) + (D5 != 0) +
77  (D6 != 0) + (D7 != 0) + (D8 != 0) &&
78  1 == D1 + D2 + D3 + D4 + D5 + D6 + D7 + D8,
79  };
80 
81  template <int R1, int R2, int R3, int R4, int R5, int R6, int R7, int R8>
82  constexpr bool operator==(dimensions<R1, R2, R3, R4, R5, R6, R7, R8> const&) const {
83  return D1 == R1 && D2 == R2 && D3 == R3 && D4 == R4 && D5 == R5 && D6 == R6 &&
84  D7 == R7 && D8 == R8;
85  }
86 
87  template <int R1, int R2, int R3, int R4, int R5, int R6, int R7, int R8>
88  constexpr bool operator!=(
90  return !(*this == rhs);
91  }
92  };
93 
95 
97 
99 
100  namespace detail {
101 
106  template <typename D, typename T>
107  struct collapse {
108  typedef quantity<D, T> type;
109  };
110 
111  template <typename T>
112  struct collapse<dimensionless_d, T> {
113  typedef T type;
114  };
115 
116  template <typename D, typename T>
117  using Collapse = typename collapse<D, T>::type;
118 
119  // promote types of expression to result type.
120 
121  template <typename X, typename Y>
122  using PromoteAdd = decltype(std::declval<X>() + std::declval<Y>());
123 
124  template <typename X, typename Y>
125  using PromoteMul = decltype(std::declval<X>() * std::declval<Y>());
126 
127  /*
128  * The following batch of structs are type generators to calculate
129  * the correct type of the result of various operations.
130  */
131 
135  template <typename DX, typename DY, typename T>
136  struct product {
137  enum {
138  d1 = DX::dim1 + DY::dim1,
139  d2 = DX::dim2 + DY::dim2,
140  d3 = DX::dim3 + DY::dim3,
141  d4 = DX::dim4 + DY::dim4,
142  d5 = DX::dim5 + DY::dim5,
143  d6 = DX::dim6 + DY::dim6,
144  d7 = DX::dim7 + DY::dim7,
145  d8 = DX::dim8 + DY::dim8
146  };
147 
149 
150  typedef Collapse<dim, T> type;
151  };
152 
157  template <typename DX, typename DY>
159 
160  template <typename DX, typename DY, typename X, typename Y>
161  using Product = typename product<DX, DY, PromoteMul<X, Y>>::type;
162 
166  template <typename DX, typename DY, typename T>
167  struct quotient {
168  enum {
169  d1 = DX::dim1 - DY::dim1,
170  d2 = DX::dim2 - DY::dim2,
171  d3 = DX::dim3 - DY::dim3,
172  d4 = DX::dim4 - DY::dim4,
173  d5 = DX::dim5 - DY::dim5,
174  d6 = DX::dim6 - DY::dim6,
175  d7 = DX::dim7 - DY::dim7,
176  d8 = DX::dim8 - DY::dim8
177  };
178 
180 
181  typedef Collapse<dim, T> type;
182  };
183 
188  template <typename DX, typename DY>
190 
191  template <typename DX, typename DY, typename X, typename Y>
192  using Quotient = typename quotient<DX, DY, PromoteMul<X, Y>>::type;
193 
197  template <typename D, typename T>
198  struct reciprocal {
199  enum {
200  d1 = -D::dim1,
201  d2 = -D::dim2,
202  d3 = -D::dim3,
203  d4 = -D::dim4,
204  d5 = -D::dim5,
205  d6 = -D::dim6,
206  d7 = -D::dim7,
207  d8 = -D::dim8
208  };
209 
211 
212  typedef Collapse<dim, T> type;
213  };
214 
219  template <typename DX, typename DY>
221 
222  template <typename D, typename X, typename Y>
223  using Reciprocal = typename reciprocal<D, PromoteMul<X, Y>>::type;
224 
228  template <typename D, int N, typename T>
229  struct power {
230  enum {
231  d1 = N * D::dim1,
232  d2 = N * D::dim2,
233  d3 = N * D::dim3,
234  d4 = N * D::dim4,
235  d5 = N * D::dim5,
236  d6 = N * D::dim6,
237  d7 = N * D::dim7,
238  d8 = N * D::dim8
239  };
240 
242 
243  typedef Collapse<dim, T> type;
244  };
245 
250  template <typename DX, int N>
251  using power_d = typename power<DX, N, Rep>::dim;
252 
253  template <typename D, int N, typename T>
254  using Power = typename power<D, N, T>::type;
255 
259  template <typename D, int N, typename T>
260  struct root {
261  enum {
262  all_even_multiples = D::dim1 % N == 0 && D::dim2 % N == 0 && D::dim3 % N == 0 &&
263  D::dim4 % N == 0 && D::dim5 % N == 0 && D::dim6 % N == 0 &&
264  D::dim7 % N == 0 && D::dim8 % N == 0
265  };
266 
267  enum {
268  d1 = D::dim1 / N,
269  d2 = D::dim2 / N,
270  d3 = D::dim3 / N,
271  d4 = D::dim4 / N,
272  d5 = D::dim5 / N,
273  d6 = D::dim6 / N,
274  d7 = D::dim7 / N,
275  d8 = D::dim8 / N
276  };
277 
279 
280  typedef Collapse<dim, T> type;
281  };
282 
286  template <typename D, int N>
287  using root_d = typename root<D, N, Rep>::dim;
288 
289  template <typename D, int N, typename T>
290  using Root = typename detail::root<D, N, T>::type;
291 
295  constexpr struct magnitude_tag_t { } magnitude_tag{}; } // namespace detail
296 
302  template <typename Dims, typename T /*= Rep */>
303  class quantity {
304  public:
305  typedef Dims dimension_type;
306 
307  typedef T value_type;
308 
309  typedef quantity<Dims, T> this_type;
310 
311  constexpr quantity()
312  : m_value{} {}
313 
318  template <typename X>
319  constexpr explicit quantity(detail::magnitude_tag_t, X x)
320  : m_value(x) {}
321 
325  template <typename X>
326  constexpr quantity(quantity<Dims, X> const& x)
327  : m_value(x.magnitude()) {}
328 
329  // /**
330  // * convert to compatible unit, for example: (3._dm).to(meter) gives 0.3;
331  // */
332  // constexpr value_type to( quantity const & x ) const { return *this / x; }
333 
337  template <typename DX, typename X>
338  constexpr auto to(quantity<DX, X> const& x) const
339  -> detail::Quotient<Dims, DX, T, X> {
340  return *this / x;
341  }
342 
346  constexpr value_type magnitude() const { return m_value; }
347 
351  constexpr dimension_type dimension() const { return dimension_type{}; }
352 
359  static constexpr quantity zero() { return quantity{value_type(0.0)}; }
360  // static constexpr quantity zero = quantity{ value_type( 0.0 ) };
361 
367  static constexpr quantity infinity() {
368  return quantity{value_type(std::numeric_limits<value_type>::infinity())};
369  }
370 
371  template <typename _dim = dimension_type,
372  std::enable_if_t<std::is_same_v<_dim, dimensionless_d>, bool> = true>
373  operator double() const {
374  return m_value;
375  }
376 
377  private:
381  constexpr explicit quantity(value_type x)
382  : m_value{x} {}
383 
384  private:
385  value_type m_value;
386 
387  enum { has_dimension = !Dims::is_all_zero };
388 
389  // static_assert( has_dimension, "quantity dimensions must not all be zero" ); //
390  // MR: removed
391 
392  private:
393  // friends:
394 
395  // arithmetic
396 
397  template <typename D, typename X, typename Y>
398  friend constexpr quantity<D, X>& operator+=(quantity<D, X>& x,
399  quantity<D, Y> const& y);
400 
401  template <typename D, typename X>
402  friend constexpr quantity<D, X> operator+(quantity<D, X> const& x);
403 
404  template <typename D, typename X, typename Y>
406  quantity<D, X> const& x, quantity<D, Y> const& y);
407 
408  template <typename D, typename X, typename Y>
409  friend constexpr quantity<D, X>& operator-=(quantity<D, X>& x,
410  quantity<D, Y> const& y);
411 
412  template <typename D, typename X>
413  friend constexpr quantity<D, X> operator-(quantity<D, X> const& x);
414 
415  template <typename D, typename X, typename Y>
417  quantity<D, X> const& x, quantity<D, Y> const& y);
418 
419  template <typename D, typename X, typename Y>
420  friend constexpr quantity<D, X>& operator*=(quantity<D, X>& x, const Y& y);
421 
422  template <typename D, typename X, typename Y>
424  quantity<D, X> const& x, const Y& y);
425 
426  template <typename D, typename X, typename Y>
428  const X& x, quantity<D, Y> const& y);
429 
430  template <typename DX, typename DY, typename X, typename Y>
431  friend constexpr detail::Product<DX, DY, X, Y> operator*(
432  quantity<DX, X> const& lhs, quantity<DY, Y> const& rhs);
433 
434  template <typename D, typename X, typename Y>
435  friend constexpr quantity<D, X>& operator/=(quantity<D, X>& x, const Y& y);
436 
437  template <typename D, typename X, typename Y>
439  quantity<D, X> const& x, const Y& y);
440 
441  template <typename D, typename X, typename Y>
442  friend constexpr detail::Reciprocal<D, X, Y> operator/(const X& x,
443  quantity<D, Y> const& y);
444 
445  template <typename DX, typename DY, typename X, typename Y>
446  friend constexpr detail::Quotient<DX, DY, X, Y> operator/(quantity<DX, X> const& x,
447  quantity<DY, Y> const& y);
448 
449  // absolute value.
450 
451  template <typename D, typename X>
452  friend constexpr quantity<D, X> abs(quantity<D, X> const& x);
453 
454  // powers and roots
455 
456  template <int N, typename D, typename X>
457  friend constexpr detail::Power<D, N, X> nth_power(quantity<D, X> const& x);
458 
459  template <typename D, typename X>
460  friend constexpr detail::Power<D, 2, X> square(quantity<D, X> const& x);
461 
462  template <typename D, typename X>
463  friend constexpr detail::Power<D, 3, X> cube(quantity<D, X> const& x);
464 
465  template <int N, typename D, typename X>
466  friend detail::Root<D, N, X> constexpr nth_root(quantity<D, X> const& x);
467 
468  template <typename D, typename X>
469  friend detail::Root<D, 2, X> constexpr sqrt(quantity<D, X> const& x);
470 
471  template <typename D, typename X>
472  friend detail::Root<D, 3, X> constexpr cbrt(quantity<D, X> const& x);
473 
474  // comparison
475 
476  template <typename D, typename X, typename Y>
477  friend constexpr bool operator==(quantity<D, X> const& x, quantity<D, Y> const& y);
478 
479  template <typename D, typename X, typename Y>
480  friend constexpr bool operator!=(quantity<D, X> const& x, quantity<D, Y> const& y);
481 
482  template <typename D, typename X, typename Y>
483  friend constexpr bool operator<(quantity<D, X> const& x, quantity<D, Y> const& y);
484 
485  template <typename D, typename X, typename Y>
486  friend constexpr bool operator<=(quantity<D, X> const& x, quantity<D, Y> const& y);
487 
488  template <typename D, typename X, typename Y>
489  friend constexpr bool operator>(quantity<D, X> const& x, quantity<D, Y> const& y);
490 
491  template <typename D, typename X, typename Y>
492  friend constexpr bool operator>=(quantity<D, X> const& x, quantity<D, Y> const& y);
493  }; // namespace units
494 
495  // helper to check whether some type is a quantity
496  template <typename TNonQuantityType>
497  struct is_quantity : std::false_type {};
498 
499  template <typename Dims, typename T>
500  struct is_quantity<quantity<Dims, T>> : std::true_type {};
501 
502  template <typename T>
503  constexpr bool is_quantity_v = is_quantity<T>::value;
504 
505  // Give names to the seven fundamental dimensions of physical reality.
506 
514  typedef dimensions<0, 0, 0, 0, 0, 0, 0, 1> hepenergy_d; // this is not an SI unit !
515 
516  // Addition operators
517 
519 
520  template <typename D, typename X, typename Y>
522  return x.m_value += y.m_value, x;
523  }
524 
526 
527  template <typename D, typename X>
529  return quantity<D, X>(+x.m_value);
530  }
531 
533 
534  template <typename D, typename X, typename Y>
536  quantity<D, Y> const& y) {
537  return quantity<D, detail::PromoteAdd<X, Y>>(x.m_value + y.m_value);
538  }
539 
540  // Subtraction operators
541 
543 
544  template <typename D, typename X, typename Y>
546  return x.m_value -= y.m_value, x;
547  }
548 
550 
551  template <typename D, typename X>
553  return quantity<D, X>(-x.m_value);
554  }
555 
557 
558  template <typename D, typename X, typename Y>
560  quantity<D, Y> const& y) {
561  return quantity<D, detail::PromoteAdd<X, Y>>(x.m_value - y.m_value);
562  }
563 
564  // Multiplication operators
565 
567 
568  template <typename D, typename X, typename Y>
569  constexpr quantity<D, X>& operator*=(quantity<D, X>& x, const Y& y) {
570  return x.m_value *= y, x;
571  }
572 
574 
575  template <typename D, typename X, typename Y>
577  const Y& y) {
578  return quantity<D, detail::PromoteMul<X, Y>>(x.m_value * y);
579  }
580 
582 
583  template <typename D, typename X, typename Y>
585  quantity<D, Y> const& y) {
586  return quantity<D, detail::PromoteMul<X, Y>>(x * y.m_value);
587  }
588 
590 
591  template <typename DX, typename DY, typename X, typename Y>
592  constexpr detail::Product<DX, DY, X, Y> operator*(quantity<DX, X> const& lhs,
593  quantity<DY, Y> const& rhs) {
594  return detail::Product<DX, DY, X, Y>(lhs.m_value * rhs.m_value);
595  }
596 
597  // Division operators
598 
600 
601  template <typename D, typename X, typename Y>
602  constexpr quantity<D, X>& operator/=(quantity<D, X>& x, const Y& y) {
603  return x.m_value /= y, x;
604  }
605 
607 
608  template <typename D, typename X, typename Y>
610  const Y& y) {
611  return quantity<D, detail::PromoteMul<X, Y>>(x.m_value / y);
612  }
613 
615 
616  template <typename D, typename X, typename Y>
617  constexpr detail::Reciprocal<D, X, Y> operator/(const X& x, quantity<D, Y> const& y) {
618  return detail::Reciprocal<D, X, Y>(x / y.m_value);
619  }
620 
622 
623  template <typename DX, typename DY, typename X, typename Y>
624  constexpr detail::Quotient<DX, DY, X, Y> operator/(quantity<DX, X> const& x,
625  quantity<DY, Y> const& y) {
626  return detail::Quotient<DX, DY, X, Y>(x.m_value / y.m_value);
627  }
628 
630 
631  template <typename D, typename X>
632  quantity<D, X> constexpr abs(quantity<D, X> const& x) {
633  return quantity<D, X>(std::abs(x.m_value));
634  }
635 
636  // General powers
637 
639 
640  template <int N, typename D, typename X>
641  detail::Power<D, N, X> constexpr nth_power(quantity<D, X> const& x) {
642  return detail::Power<D, N, X>(std::pow(x.m_value, X(N)));
643  }
644 
645  // Low powers defined separately for efficiency.
646 
648 
649  template <typename D, typename X>
650  constexpr detail::Power<D, 2, X> square(quantity<D, X> const& x) {
651  return x * x;
652  }
653 
655 
656  template <typename D, typename X>
657  constexpr detail::Power<D, 3, X> cube(quantity<D, X> const& x) {
658  return x * x * x;
659  }
660 
661  // General root
662 
664 
665  template <int N, typename D, typename X>
666  detail::Root<D, N, X> constexpr nth_root(quantity<D, X> const& x) {
668  "root result dimensions must be integral");
669 
670  return detail::Root<D, N, X>(std::pow(x.m_value, X(1.0) / N));
671  }
672 
673  // Low roots defined separately for convenience.
674 
676 
677  template <typename D, typename X>
678  detail::Root<D, 2, X> constexpr sqrt(quantity<D, X> const& x) {
680  "root result dimensions must be integral");
681 
682  return detail::Root<D, 2, X>(std::sqrt(x.m_value));
683  }
684 
686 
687  template <typename D, typename X>
688  detail::Root<D, 3, X> constexpr cbrt(quantity<D, X> const& x) {
690  "root result dimensions must be integral");
691 
692  return detail::Root<D, 3, X>(std::cbrt(x.m_value));
693  }
694 
695  // Comparison operators
696 
698 
699  template <typename D, typename X, typename Y>
700  constexpr bool operator==(quantity<D, X> const& x, quantity<D, Y> const& y) {
701  return x.m_value == y.m_value;
702  }
703 
705 
706  template <typename D, typename X, typename Y>
707  constexpr bool operator!=(quantity<D, X> const& x, quantity<D, Y> const& y) {
708  return x.m_value != y.m_value;
709  }
710 
712 
713  template <typename D, typename X, typename Y>
714  constexpr bool operator<(quantity<D, X> const& x, quantity<D, Y> const& y) {
715  return x.m_value < y.m_value;
716  }
717 
719 
720  template <typename D, typename X, typename Y>
721  constexpr bool operator<=(quantity<D, X> const& x, quantity<D, Y> const& y) {
722  return x.m_value <= y.m_value;
723  }
724 
726 
727  template <typename D, typename X, typename Y>
728  constexpr bool operator>(quantity<D, X> const& x, quantity<D, Y> const& y) {
729  return x.m_value > y.m_value;
730  }
731 
733 
734  template <typename D, typename X, typename Y>
735  constexpr bool operator>=(quantity<D, X> const& x, quantity<D, Y> const& y) {
736  return x.m_value >= y.m_value;
737  }
738 
740 
741  template <typename DX, typename X>
742  inline constexpr DX dimension(quantity<DX, X> const& q) {
743  return q.dimension();
744  }
745 
747 
748  template <typename DX, typename X>
749  inline constexpr X magnitude(quantity<DX, X> const& q) {
750  return q.magnitude();
751  }
752 
753  // The seven SI base units. These tie our numbers to the real world.
754 
755  constexpr quantity<length_d> meter{detail::magnitude_tag, 1.0};
756  constexpr quantity<mass_d> kilogram{detail::magnitude_tag, 1.0};
757  constexpr quantity<time_interval_d> second{detail::magnitude_tag, 1.0};
758  constexpr quantity<electric_current_d> ampere{detail::magnitude_tag, 1.0};
759  constexpr quantity<thermodynamic_temperature_d> kelvin{detail::magnitude_tag, 1.0};
760  constexpr quantity<amount_of_substance_d> mole{detail::magnitude_tag, 1.0};
761  constexpr quantity<luminous_intensity_d> candela{detail::magnitude_tag, 1.0};
762  constexpr quantity<hepenergy_d> electronvolt{detail::magnitude_tag, 1.0};
763 
764  // The standard SI prefixes.
765 
766  constexpr long double quetta = 1e+30L;
767  constexpr long double ronna = 1e+27L;
768  constexpr long double yotta = 1e+24L;
769  constexpr long double zetta = 1e+21L;
770  constexpr long double exa = 1e+18L;
771  constexpr long double peta = 1e+15L;
772  constexpr long double tera = 1e+12L;
773  constexpr long double giga = 1e+9L;
774  constexpr long double mega = 1e+6L;
775  constexpr long double kilo = 1e+3L;
776  constexpr long double hecto = 1e+2L;
777  constexpr long double deka = 1e+1L;
778  constexpr long double deci = 1e-1L;
779  constexpr long double centi = 1e-2L;
780  constexpr long double milli = 1e-3L;
781  constexpr long double micro = 1e-6L;
782  constexpr long double nano = 1e-9L;
783  constexpr long double pico = 1e-12L;
784  constexpr long double femto = 1e-15L;
785  constexpr long double atto = 1e-18L;
786  constexpr long double zepto = 1e-21L;
787  constexpr long double yocto = 1e-24L;
788  constexpr long double ronto = 1e-27L;
789  constexpr long double quecto = 1e-30L;
790 
791  // Binary prefixes, pending adoption.
792 
793  constexpr long double kibi = 1024;
794  constexpr long double mebi = 1024 * kibi;
795  constexpr long double gibi = 1024 * mebi;
796  constexpr long double tebi = 1024 * gibi;
797  constexpr long double pebi = 1024 * tebi;
798  constexpr long double exbi = 1024 * pebi;
799  constexpr long double zebi = 1024 * exbi;
800  constexpr long double yobi = 1024 * zebi;
801 
802  // The rest of the standard dimensional types, as specified in SP811.
803 
804  using absorbed_dose_d = dimensions<2, 0, -2>;
805  using absorbed_dose_rate_d = dimensions<2, 0, -3>;
806  using acceleration_d = dimensions<1, 0, -2>;
807  using activity_of_a_nuclide_d = dimensions<0, 0, -1>;
808  using angular_velocity_d = dimensions<0, 0, -1>;
809  using angular_acceleration_d = dimensions<0, 0, -2>;
810  using area_d = dimensions<2, 0, 0>;
811  using capacitance_d = dimensions<-2, -1, 4, 2>;
812  using concentration_d = dimensions<-3, 0, 0, 0, 0, 1>;
813  using current_density_d = dimensions<-2, 0, 0, 1>;
814  using dose_equivalent_d = dimensions<2, 0, -2>;
815  using dynamic_viscosity_d = dimensions<-1, 1, -1>;
817  using electric_charge_density_d = dimensions<-3, 0, 1, 1>;
818  using electric_conductance_d = dimensions<-2, -1, 3, 2>;
819  using electric_field_strenth_d = dimensions<1, 1, -3, -1>;
820  using electric_flux_density_d = dimensions<-2, 0, 1, 1>;
821  using electric_potential_d = dimensions<2, 1, -3, -1>;
822  using electric_resistance_d = dimensions<2, 1, -3, -2>;
823  using energy_d = dimensions<2, 1, -2>;
824  using energy_density_d = dimensions<-1, 1, -2>;
825  using exposure_d = dimensions<0, -1, 1, 1>;
826  using force_d = dimensions<1, 1, -2>;
827  using frequency_d = dimensions<0, 0, -1>;
828  using heat_capacity_d = dimensions<2, 1, -2, 0, -1>;
829  using heat_density_d = dimensions<0, 1, -2>;
830  using heat_density_flow_rate_d = dimensions<0, 1, -3>;
831  using heat_flow_rate_d = dimensions<2, 1, -3>;
832  using heat_flux_density_d = dimensions<0, 1, -3>;
833  using heat_transfer_coefficient_d = dimensions<0, 1, -3, 0, -1>;
834  using illuminance_d = dimensions<-2, 0, 0, 0, 0, 0, 1>;
835  using inductance_d = dimensions<2, 1, -2, -2>;
836  using irradiance_d = dimensions<0, 1, -3>;
837  using kinematic_viscosity_d = dimensions<2, 0, -1>;
838  using luminance_d = dimensions<-2, 0, 0, 0, 0, 0, 1>;
840  using magnetic_field_strength_d = dimensions<-1, 0, 0, 1>;
841  using magnetic_flux_d = dimensions<2, 1, -2, -1>;
842  using magnetic_flux_density_d = dimensions<0, 1, -2, -1>;
843  using magnetic_permeability_d = dimensions<1, 1, -2, -2>;
844  using mass_density_d = dimensions<-3, 1, 0>;
845  using mass_flow_rate_d = dimensions<0, 1, -1>;
846  using molar_energy_d = dimensions<2, 1, -2, 0, 0, -1>;
847  using molar_entropy_d = dimensions<2, 1, -2, -1, 0, -1>;
848  using moment_of_force_d = dimensions<2, 1, -2>;
849  using permittivity_d = dimensions<-3, -1, 4, 2>;
850  using power_d = dimensions<2, 1, -3>;
851  using pressure_d = dimensions<-1, 1, -2>;
852  using radiance_d = dimensions<0, 1, -3>;
853  using radiant_intensity_d = dimensions<2, 1, -3>;
854  using speed_d = dimensions<1, 0, -1>;
855  using specific_energy_d = dimensions<2, 0, -2>;
856  using specific_heat_capacity_d = dimensions<2, 0, -2, 0, -1>;
857  using specific_volume_d = dimensions<3, -1, 0>;
858  using substance_permeability_d = dimensions<-1, 0, 1>;
859  using surface_tension_d = dimensions<0, 1, -2>;
860  using thermal_conductivity_d = dimensions<1, 1, -3, 0, -1>;
861  using thermal_diffusivity_d = dimensions<2, 0, -1>;
862  using thermal_insulance_d = dimensions<0, -1, 3, 0, 1>;
863  using thermal_resistance_d = dimensions<-2, -1, 3, 0, 1>;
864  using thermal_resistivity_d = dimensions<-1, -1, 3, 0, 1>;
865  using torque_d = dimensions<2, 1, -2>;
867  using volume_flow_rate_d = dimensions<3, 0, -1>;
868  using wave_number_d = dimensions<-1, 0, 0>;
869 
870  // Handy values.
871 
872  constexpr Rep pi{Rep(3.141592653589793238462L)};
873  constexpr Rep percent{Rep(1) / 100};
874 
876  constexpr quantity<mass_d> gram{kilogram / 1000};
877 
878  // The derived SI units, as specified in SP811.
879 
880  constexpr Rep radian{Rep(1)};
881  constexpr Rep steradian{Rep(1)};
882  constexpr quantity<force_d> newton{meter * kilogram / square(second)};
883  constexpr quantity<pressure_d> pascal{newton / square(meter)};
884  constexpr quantity<energy_d> joule{newton * meter};
885  constexpr quantity<power_d> watt{joule / second};
886  constexpr quantity<electric_charge_d> coulomb{second * ampere};
887  constexpr quantity<electric_potential_d> volt{watt / ampere};
888  constexpr quantity<capacitance_d> farad{coulomb / volt};
889  constexpr quantity<electric_resistance_d> ohm{volt / ampere};
890  constexpr quantity<electric_conductance_d> siemens{ampere / volt};
891  constexpr quantity<magnetic_flux_d> weber{volt * second};
892  constexpr quantity<magnetic_flux_density_d> tesla{weber / square(meter)};
893  constexpr quantity<inductance_d> henry{weber / ampere};
894  constexpr quantity<thermodynamic_temperature_d> degree_celsius{kelvin};
895  constexpr quantity<luminous_flux_d> lumen{candela * steradian};
896  constexpr quantity<illuminance_d> lux{lumen / meter / meter};
897  constexpr quantity<activity_of_a_nuclide_d> becquerel{1 / second};
898  constexpr quantity<absorbed_dose_d> gray{joule / kilogram};
899  constexpr quantity<dose_equivalent_d> sievert{joule / kilogram};
900  constexpr quantity<frequency_d> hertz{1 / second};
901 
902  // The rest of the units approved for use with SI, as specified in SP811.
903  // (However, use of these units is generally discouraged.)
904 
905  constexpr quantity<length_d> angstrom{Rep(1e-10L) * meter};
906  constexpr quantity<area_d> are{Rep(1e+2L) * square(meter)};
907  constexpr quantity<pressure_d> bar{Rep(1e+5L) * pascal};
908  constexpr quantity<area_d> barn{Rep(1e-28L) * square(meter)};
909  constexpr quantity<activity_of_a_nuclide_d> curie{Rep(3.7e+10L) * becquerel};
910  constexpr quantity<time_interval_d> day{Rep(86400L) * second};
911  constexpr Rep degree_angle{pi / 180};
912  constexpr quantity<acceleration_d> gal{Rep(1e-2L) * meter / square(second)};
913  constexpr quantity<area_d> hectare{Rep(1e+4L) * square(meter)};
914  constexpr quantity<time_interval_d> hour{Rep(3600) * second};
915  constexpr quantity<speed_d> knot{Rep(1852) / 3600 * meter / second};
916  constexpr quantity<volume_d> liter{Rep(1e-3L) * cube(meter)};
917  constexpr quantity<time_interval_d> minute{Rep(60) * second};
918  constexpr Rep minute_angle{pi / 10800};
919  constexpr quantity<length_d> mile_nautical{Rep(1852) * meter};
920  constexpr quantity<absorbed_dose_d> rad{Rep(1e-2L) * gray};
921  constexpr quantity<dose_equivalent_d> rem{Rep(1e-2L) * sievert};
922  constexpr quantity<exposure_d> roentgen{Rep(2.58e-4L) * coulomb / kilogram};
923  constexpr Rep second_angle{pi / 648000L};
924  constexpr quantity<mass_d> ton_metric{Rep(1e+3L) * kilogram};
925 
926  // Alternate (non-US) spellings:
927 
928  constexpr quantity<length_d> metre{meter};
929  constexpr quantity<volume_d> litre{liter};
930  constexpr Rep deca{deka};
931  constexpr quantity<mass_d> tonne{ton_metric};
932 
933  // cooked literals for base units;
934  // these could also have been created with a script.
935 
936 #define QUANTITY_DEFINE_SCALING_LITERAL(sfx, dim, factor) \
937  constexpr quantity<dim, double> operator"" _##sfx(unsigned long long x) { \
938  return quantity<dim, double>(detail::magnitude_tag, factor * x); \
939  } \
940  constexpr quantity<dim, double> operator"" _##sfx(long double x) { \
941  return quantity<dim, double>(detail::magnitude_tag, factor * x); \
942  }
943 
944 #define QUANTITY_DEFINE_SCALING_LITERALS(pfx, dim, fact) \
945  QUANTITY_DEFINE_SCALING_LITERAL(Q##pfx, dim, fact* quetta) \
946  QUANTITY_DEFINE_SCALING_LITERAL(R##pfx, dim, fact* ronna) \
947  QUANTITY_DEFINE_SCALING_LITERAL(Y##pfx, dim, fact* yotta) \
948  QUANTITY_DEFINE_SCALING_LITERAL(Z##pfx, dim, fact* zetta) \
949  QUANTITY_DEFINE_SCALING_LITERAL(E##pfx, dim, fact* exa) \
950  QUANTITY_DEFINE_SCALING_LITERAL(P##pfx, dim, fact* peta) \
951  QUANTITY_DEFINE_SCALING_LITERAL(T##pfx, dim, fact* tera) \
952  QUANTITY_DEFINE_SCALING_LITERAL(G##pfx, dim, fact* giga) \
953  QUANTITY_DEFINE_SCALING_LITERAL(M##pfx, dim, fact* mega) \
954  QUANTITY_DEFINE_SCALING_LITERAL(k##pfx, dim, fact* kilo) \
955  QUANTITY_DEFINE_SCALING_LITERAL(h##pfx, dim, fact* hecto) \
956  QUANTITY_DEFINE_SCALING_LITERAL(da##pfx, dim, fact* deka) \
957  QUANTITY_DEFINE_SCALING_LITERAL(pfx, dim, fact * 1) \
958  QUANTITY_DEFINE_SCALING_LITERAL(d##pfx, dim, fact* deci) \
959  QUANTITY_DEFINE_SCALING_LITERAL(c##pfx, dim, fact* centi) \
960  QUANTITY_DEFINE_SCALING_LITERAL(m##pfx, dim, fact* milli) \
961  QUANTITY_DEFINE_SCALING_LITERAL(u##pfx, dim, fact* micro) \
962  QUANTITY_DEFINE_SCALING_LITERAL(n##pfx, dim, fact* nano) \
963  QUANTITY_DEFINE_SCALING_LITERAL(p##pfx, dim, fact* pico) \
964  QUANTITY_DEFINE_SCALING_LITERAL(f##pfx, dim, fact* femto) \
965  QUANTITY_DEFINE_SCALING_LITERAL(a##pfx, dim, fact* atto) \
966  QUANTITY_DEFINE_SCALING_LITERAL(z##pfx, dim, fact* zepto) \
967  QUANTITY_DEFINE_SCALING_LITERAL(y##pfx, dim, fact* yocto) \
968  QUANTITY_DEFINE_SCALING_LITERAL(r##pfx, dim, fact* ronto) \
969  QUANTITY_DEFINE_SCALING_LITERAL(q##pfx, dim, fact* quecto)
970 
971 #define QUANTITY_DEFINE_LITERALS(pfx, dim) QUANTITY_DEFINE_SCALING_LITERALS(pfx, dim, 1)
972 
974 
975  namespace literals {
976 
977  QUANTITY_DEFINE_SCALING_LITERALS(g, mass_d, 1e-3)
978 
979  QUANTITY_DEFINE_LITERALS(m, length_d)
980  QUANTITY_DEFINE_LITERALS(s, time_interval_d)
981  QUANTITY_DEFINE_LITERALS(A, electric_current_d)
982  QUANTITY_DEFINE_LITERALS(K, thermodynamic_temperature_d)
983  QUANTITY_DEFINE_LITERALS(mol, amount_of_substance_d)
984  QUANTITY_DEFINE_LITERALS(cd, luminous_intensity_d)
985 
986  } // namespace literals
987 
988  } // namespace units
989 } // namespace phys
990 
991 #endif // PHYS_UNITS_QUANTITY_HPP_INCLUDED
992 
993 /*
994  * end of file
995  */
quantity< D, X > constexpr abs(quantity< D, X > const &x)
absolute value.
Definition: quantity.hpp:632
We could drag dimensions around individually, but it&#39;s much more convenient to package them...
Definition: quantity.hpp:62
detail::Root< D, N, X > constexpr nth_root(quantity< D, X > const &x)
n-th root.
Definition: quantity.hpp:666
constexpr auto to(quantity< DX, X > const &x) const -> detail::Quotient< Dims, DX, T, X >
convert to compatible unit, for example: (3._dm).to(meter) gives 0.3;
Definition: quantity.hpp:338
constexpr quantity< D, X > operator-(quantity< D, X > const &x)
Definition: quantity.hpp:552
dimensions< 0, 0, 0 > dimensionless_d
demensionless &#39;dimension&#39;.
Definition: quantity.hpp:96
product type generator.
Definition: quantity.hpp:136
constexpr quantity< D, detail::PromoteMul< X, Y > > operator*(quantity< D, X > const &x, const Y &y)
quan * num
Definition: quantity.hpp:576
namespace phys.
constexpr quantity(detail::magnitude_tag_t, X x)
public converting initializing constructor; requires magnitude_tag to prevent constructing a quantity...
Definition: quantity.hpp:319
reciprocal type generator.
Definition: quantity.hpp:198
constexpr quantity< D, detail::PromoteMul< X, Y > > operator/(quantity< D, X > const &x, const Y &y)
quan / num
Definition: quantity.hpp:609
The "collapse" template is used to avoid quantity< dimensions< 0, 0, 0 > >, i.e.
Definition: quantity.hpp:107
constexpr bool operator>=(quantity< D, X > const &x, quantity< D, Y > const &y)
greater-equal.
Definition: quantity.hpp:735
constexpr dimension_type dimension() const
the quantity&#39;s dimensions.
Definition: quantity.hpp:351
detail::Root< D, 3, X > constexpr cbrt(quantity< D, X > const &x)
cubic root.
Definition: quantity.hpp:688
constexpr detail::Power< D, 3, X > cube(quantity< D, X > const &x)
cube.
Definition: quantity.hpp:657
constexpr quantity< D, X > operator+(quantity< D, X > const &x)
Definition: quantity.hpp:528
static constexpr quantity infinity()
We also define "infinity" for each type.
Definition: quantity.hpp:367
constexpr quantity(quantity< Dims, X > const &x)
converting copy-assignment constructor.
Definition: quantity.hpp:326
typename reciprocal< DX, Rep >::dim reciprocal_d
convenience type for reciprocal dimension
Definition: quantity.hpp:220
constexpr DX dimension(quantity< DX, X > const &q)
quantity&#39;s dimension.
Definition: quantity.hpp:742
class "quantity" is the heart of the library.
Definition: quantity.hpp:54
tag to construct a quantity from a magnitude.
Definition: quantity.hpp:295
constexpr value_type magnitude() const
the quantity&#39;s magnitude.
Definition: quantity.hpp:346
constexpr detail::Power< D, 2, X > square(quantity< D, X > const &x)
square.
Definition: quantity.hpp:650
static constexpr quantity zero()
We need a "zero" of each type – for comparisons, to initialize running totals, etc.
Definition: quantity.hpp:359
typename root< D, N, Rep >::dim root_d
convenience type for root dimension
Definition: quantity.hpp:287
typename product< DX, DY, Rep >::dim product_d
convenience type for product dimension
Definition: quantity.hpp:158
constexpr quantity< D, X > & operator-=(quantity< D, X > &x, quantity< D, Y > const &y)
quan -= quan
Definition: quantity.hpp:545
detail::Root< D, 2, X > constexpr sqrt(quantity< D, X > const &x)
square root.
Definition: quantity.hpp:678
typename power< DX, N, Rep >::dim power_d
convenience type for power dimension
Definition: quantity.hpp:251
power type generator.
Definition: quantity.hpp:229
constexpr X magnitude(quantity< DX, X > const &q)
quantity&#39;s magnitude.
Definition: quantity.hpp:749
quotient type generator.
Definition: quantity.hpp:167
constexpr quantity< D, X > & operator*=(quantity< D, X > &x, const Y &y)
quan *= num
Definition: quantity.hpp:569
root type generator.
Definition: quantity.hpp:260
detail::Power< D, N, X > constexpr nth_power(quantity< D, X > const &x)
N-th power.
Definition: quantity.hpp:641
constexpr quantity< D, X > & operator+=(quantity< D, X > &x, quantity< D, Y > const &y)
quan += quan
Definition: quantity.hpp:521
typename quotient< DX, DY, Rep >::dim quotient_d
convenience type for quotient dimension
Definition: quantity.hpp:189
constexpr bool operator>(quantity< D, X > const &x, quantity< D, Y > const &y)
greater-than.
Definition: quantity.hpp:728
constexpr quantity< D, X > & operator/=(quantity< D, X > &x, const Y &y)
quan /= num
Definition: quantity.hpp:602