CORSIKA  @c8_version@
The framework to simulate particle cascades for astroparticle physics
quantity.hpp File Reference

Zero-overhead dimensional analysis and unit/quantity manipulation and conversion. More...

#include <cmath>
#include <cstdlib>
#include <utility>
#include <type_traits>
Include dependency graph for quantity.hpp:
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Classes

class  phys::units::quantity< Dims, T >
 class "quantity" is the heart of the library. More...
 
struct  phys::units::dimensions< D1, D2, D3, D4, D5, D6, D7, D8 >
 We could drag dimensions around individually, but it's much more convenient to package them. More...
 
struct  phys::units::detail::collapse< D, T >
 The "collapse" template is used to avoid quantity< dimensions< 0, 0, 0 > >, i.e. More...
 
struct  phys::units::detail::collapse< dimensionless_d, T >
 
struct  phys::units::detail::product< DX, DY, T >
 product type generator. More...
 
struct  phys::units::detail::quotient< DX, DY, T >
 quotient type generator. More...
 
struct  phys::units::detail::reciprocal< D, T >
 reciprocal type generator. More...
 
struct  phys::units::detail::power< D, N, T >
 power type generator. More...
 
struct  phys::units::detail::root< D, N, T >
 root type generator. More...
 
struct  phys::units::detail::magnitude_tag_t
 tag to construct a quantity from a magnitude. More...
 
class  phys::units::quantity< Dims, T >
 class "quantity" is the heart of the library. More...
 
struct  phys::units::is_quantity< TNonQuantityType >
 
struct  phys::units::is_quantity< quantity< Dims, T > >
 

Namespaces

 phys
 namespace phys.
 
 phys::units
 namespace units.
 
 phys::units::detail
 namespace detail.
 
 phys::units::literals
 literals
 

Macros

#define QUANTITY_DEFINE_SCALING_LITERAL(sfx, dim, factor)
 
#define QUANTITY_DEFINE_SCALING_LITERALS(pfx, dim, fact)
 
#define QUANTITY_DEFINE_LITERALS(pfx, dim)   QUANTITY_DEFINE_SCALING_LITERALS(pfx, dim, 1)
 

Typedefs

using phys::units::Rep = double
 
typedef dimensions< 0, 0, 0 > phys::units::dimensionless_d
 demensionless 'dimension'.
 
template<typename D , typename T >
using phys::units::detail::Collapse = typename collapse< D, T >::type
 
template<typename X , typename Y >
using phys::units::detail::PromoteAdd = decltype(std::declval< X >()+std::declval< Y >())
 
template<typename X , typename Y >
using phys::units::detail::PromoteMul = decltype(std::declval< X >() *std::declval< Y >())
 
template<typename DX , typename DY >
using phys::units::detail::product_d = typename product< DX, DY, Rep >::dim
 convenience type for product dimension
 
template<typename DX , typename DY , typename X , typename Y >
using phys::units::detail::Product = typename product< DX, DY, PromoteMul< X, Y > >::type
 
template<typename DX , typename DY >
using phys::units::detail::quotient_d = typename quotient< DX, DY, Rep >::dim
 convenience type for quotient dimension
 
template<typename DX , typename DY , typename X , typename Y >
using phys::units::detail::Quotient = typename quotient< DX, DY, PromoteMul< X, Y > >::type
 
template<typename DX , typename DY >
using phys::units::detail::reciprocal_d = typename reciprocal< DX, Rep >::dim
 convenience type for reciprocal dimension
 
template<typename D , typename X , typename Y >
using phys::units::detail::Reciprocal = typename reciprocal< D, PromoteMul< X, Y > >::type
 
template<typename DX , int N>
using phys::units::detail::power_d = typename power< DX, N, Rep >::dim
 convenience type for power dimension
 
template<typename D , int N, typename T >
using phys::units::detail::Power = typename power< D, N, T >::type
 
template<typename D , int N>
using phys::units::detail::root_d = typename root< D, N, Rep >::dim
 convenience type for root dimension
 
template<typename D , int N, typename T >
using phys::units::detail::Root = typename detail::root< D, N, T >::type
 
typedef dimensions< 1, 0, 0, 0, 0, 0, 0, 0 > phys::units::length_d
 
typedef dimensions< 0, 1, 0, 0, 0, 0, 0, 0 > phys::units::mass_d
 
typedef dimensions< 0, 0, 1, 0, 0, 0, 0, 0 > phys::units::time_interval_d
 
typedef dimensions< 0, 0, 0, 1, 0, 0, 0, 0 > phys::units::electric_current_d
 
typedef dimensions< 0, 0, 0, 0, 1, 0, 0, 0 > phys::units::thermodynamic_temperature_d
 
typedef dimensions< 0, 0, 0, 0, 0, 1, 0, 0 > phys::units::amount_of_substance_d
 
typedef dimensions< 0, 0, 0, 0, 0, 0, 1, 0 > phys::units::luminous_intensity_d
 
typedef dimensions< 0, 0, 0, 0, 0, 0, 0, 1 > phys::units::hepenergy_d
 
using phys::units::absorbed_dose_d = dimensions< 2, 0, -2 >
 
using phys::units::absorbed_dose_rate_d = dimensions< 2, 0, -3 >
 
using phys::units::acceleration_d = dimensions< 1, 0, -2 >
 
using phys::units::activity_of_a_nuclide_d = dimensions< 0, 0, -1 >
 
using phys::units::angular_velocity_d = dimensions< 0, 0, -1 >
 
using phys::units::angular_acceleration_d = dimensions< 0, 0, -2 >
 
using phys::units::area_d = dimensions< 2, 0, 0 >
 
using phys::units::capacitance_d = dimensions<-2, -1, 4, 2 >
 
using phys::units::concentration_d = dimensions<-3, 0, 0, 0, 0, 1 >
 
using phys::units::current_density_d = dimensions<-2, 0, 0, 1 >
 
using phys::units::dose_equivalent_d = dimensions< 2, 0, -2 >
 
using phys::units::dynamic_viscosity_d = dimensions<-1, 1, -1 >
 
using phys::units::electric_charge_d = dimensions< 0, 0, 1, 1 >
 
using phys::units::electric_charge_density_d = dimensions<-3, 0, 1, 1 >
 
using phys::units::electric_conductance_d = dimensions<-2, -1, 3, 2 >
 
using phys::units::electric_field_strenth_d = dimensions< 1, 1, -3, -1 >
 
using phys::units::electric_flux_density_d = dimensions<-2, 0, 1, 1 >
 
using phys::units::electric_potential_d = dimensions< 2, 1, -3, -1 >
 
using phys::units::electric_resistance_d = dimensions< 2, 1, -3, -2 >
 
using phys::units::energy_d = dimensions< 2, 1, -2 >
 
using phys::units::energy_density_d = dimensions<-1, 1, -2 >
 
using phys::units::exposure_d = dimensions< 0, -1, 1, 1 >
 
using phys::units::force_d = dimensions< 1, 1, -2 >
 
using phys::units::frequency_d = dimensions< 0, 0, -1 >
 
using phys::units::heat_capacity_d = dimensions< 2, 1, -2, 0, -1 >
 
using phys::units::heat_density_d = dimensions< 0, 1, -2 >
 
using phys::units::heat_density_flow_rate_d = dimensions< 0, 1, -3 >
 
using phys::units::heat_flow_rate_d = dimensions< 2, 1, -3 >
 
using phys::units::heat_flux_density_d = dimensions< 0, 1, -3 >
 
using phys::units::heat_transfer_coefficient_d = dimensions< 0, 1, -3, 0, -1 >
 
using phys::units::illuminance_d = dimensions<-2, 0, 0, 0, 0, 0, 1 >
 
using phys::units::inductance_d = dimensions< 2, 1, -2, -2 >
 
using phys::units::irradiance_d = dimensions< 0, 1, -3 >
 
using phys::units::kinematic_viscosity_d = dimensions< 2, 0, -1 >
 
using phys::units::luminance_d = dimensions<-2, 0, 0, 0, 0, 0, 1 >
 
using phys::units::luminous_flux_d = dimensions< 0, 0, 0, 0, 0, 0, 1 >
 
using phys::units::magnetic_field_strength_d = dimensions<-1, 0, 0, 1 >
 
using phys::units::magnetic_flux_d = dimensions< 2, 1, -2, -1 >
 
using phys::units::magnetic_flux_density_d = dimensions< 0, 1, -2, -1 >
 
using phys::units::magnetic_permeability_d = dimensions< 1, 1, -2, -2 >
 
using phys::units::mass_density_d = dimensions<-3, 1, 0 >
 
using phys::units::mass_flow_rate_d = dimensions< 0, 1, -1 >
 
using phys::units::molar_energy_d = dimensions< 2, 1, -2, 0, 0, -1 >
 
using phys::units::molar_entropy_d = dimensions< 2, 1, -2, -1, 0, -1 >
 
using phys::units::moment_of_force_d = dimensions< 2, 1, -2 >
 
using phys::units::permittivity_d = dimensions<-3, -1, 4, 2 >
 
using phys::units::power_d = dimensions< 2, 1, -3 >
 
using phys::units::pressure_d = dimensions<-1, 1, -2 >
 
using phys::units::radiance_d = dimensions< 0, 1, -3 >
 
using phys::units::radiant_intensity_d = dimensions< 2, 1, -3 >
 
using phys::units::speed_d = dimensions< 1, 0, -1 >
 
using phys::units::specific_energy_d = dimensions< 2, 0, -2 >
 
using phys::units::specific_heat_capacity_d = dimensions< 2, 0, -2, 0, -1 >
 
using phys::units::specific_volume_d = dimensions< 3, -1, 0 >
 
using phys::units::substance_permeability_d = dimensions<-1, 0, 1 >
 
using phys::units::surface_tension_d = dimensions< 0, 1, -2 >
 
using phys::units::thermal_conductivity_d = dimensions< 1, 1, -3, 0, -1 >
 
using phys::units::thermal_diffusivity_d = dimensions< 2, 0, -1 >
 
using phys::units::thermal_insulance_d = dimensions< 0, -1, 3, 0, 1 >
 
using phys::units::thermal_resistance_d = dimensions<-2, -1, 3, 0, 1 >
 
using phys::units::thermal_resistivity_d = dimensions<-1, -1, 3, 0, 1 >
 
using phys::units::torque_d = dimensions< 2, 1, -2 >
 
using phys::units::volume_d = dimensions< 3, 0, 0 >
 
using phys::units::volume_flow_rate_d = dimensions< 3, 0, -1 >
 
using phys::units::wave_number_d = dimensions<-1, 0, 0 >
 

Functions

template<typename D , typename X , typename Y >
constexpr quantity< D, X > & phys::units::operator+= (quantity< D, X > &x, quantity< D, Y > const &y)
 quan += quan
 
template<typename D , typename X >
constexpr quantity< D, X > phys::units::operator+ (quantity< D, X > const &x)
 
template<typename D , typename X , typename Y >
constexpr quantity< D, detail::PromoteAdd< X, Y > > phys::units::operator+ (quantity< D, X > const &x, quantity< D, Y > const &y)
 quan + quan
 
template<typename D , typename X , typename Y >
constexpr quantity< D, X > & phys::units::operator-= (quantity< D, X > &x, quantity< D, Y > const &y)
 quan -= quan
 
template<typename D , typename X >
constexpr quantity< D, X > phys::units::operator- (quantity< D, X > const &x)
 
template<typename D , typename X , typename Y >
constexpr quantity< D, detail::PromoteAdd< X, Y > > phys::units::operator- (quantity< D, X > const &x, quantity< D, Y > const &y)
 quan - quan
 
template<typename D , typename X , typename Y >
constexpr quantity< D, X > & phys::units::operator*= (quantity< D, X > &x, const Y &y)
 quan *= num
 
template<typename D , typename X , typename Y >
constexpr quantity< D, detail::PromoteMul< X, Y > > phys::units::operator* (quantity< D, X > const &x, const Y &y)
 quan * num
 
template<typename D , typename X , typename Y >
constexpr quantity< D, detail::PromoteMul< X, Y > > phys::units::operator* (const X &x, quantity< D, Y > const &y)
 num * quan
 
template<typename DX , typename DY , typename X , typename Y >
constexpr detail::Product< DX, DY, X, Y > phys::units::operator* (quantity< DX, X > const &lhs, quantity< DY, Y > const &rhs)
 quan * quan:
 
template<typename D , typename X , typename Y >
constexpr quantity< D, X > & phys::units::operator/= (quantity< D, X > &x, const Y &y)
 quan /= num
 
template<typename D , typename X , typename Y >
constexpr quantity< D, detail::PromoteMul< X, Y > > phys::units::operator/ (quantity< D, X > const &x, const Y &y)
 quan / num
 
template<typename D , typename X , typename Y >
constexpr detail::Reciprocal< D, X, Y > phys::units::operator/ (const X &x, quantity< D, Y > const &y)
 num / quan
 
template<typename DX , typename DY , typename X , typename Y >
constexpr detail::Quotient< DX, DY, X, Y > phys::units::operator/ (quantity< DX, X > const &x, quantity< DY, Y > const &y)
 quan / quan:
 
template<typename D , typename X >
quantity< D, X > constexpr phys::units::abs (quantity< D, X > const &x)
 absolute value.
 
template<int N, typename D , typename X >
detail::Power< D, N, X > constexpr phys::units::nth_power (quantity< D, X > const &x)
 N-th power.
 
template<typename D , typename X >
constexpr detail::Power< D, 2, X > phys::units::square (quantity< D, X > const &x)
 square.
 
template<typename D , typename X >
constexpr detail::Power< D, 3, X > phys::units::cube (quantity< D, X > const &x)
 cube.
 
template<int N, typename D , typename X >
detail::Root< D, N, X > constexpr phys::units::nth_root (quantity< D, X > const &x)
 n-th root.
 
template<typename D , typename X >
detail::Root< D, 2, X > constexpr phys::units::sqrt (quantity< D, X > const &x)
 square root.
 
template<typename D , typename X >
detail::Root< D, 3, X > constexpr phys::units::cbrt (quantity< D, X > const &x)
 cubic root.
 
template<typename D , typename X , typename Y >
constexpr bool phys::units::operator== (quantity< D, X > const &x, quantity< D, Y > const &y)
 equality.
 
template<typename D , typename X , typename Y >
constexpr bool phys::units::operator!= (quantity< D, X > const &x, quantity< D, Y > const &y)
 inequality.
 
template<typename D , typename X , typename Y >
constexpr bool phys::units::operator< (quantity< D, X > const &x, quantity< D, Y > const &y)
 less-than.
 
template<typename D , typename X , typename Y >
constexpr bool phys::units::operator<= (quantity< D, X > const &x, quantity< D, Y > const &y)
 less-equal.
 
template<typename D , typename X , typename Y >
constexpr bool phys::units::operator> (quantity< D, X > const &x, quantity< D, Y > const &y)
 greater-than.
 
template<typename D , typename X , typename Y >
constexpr bool phys::units::operator>= (quantity< D, X > const &x, quantity< D, Y > const &y)
 greater-equal.
 
template<typename DX , typename X >
constexpr DX phys::units::dimension (quantity< DX, X > const &q)
 quantity's dimension.
 
template<typename DX , typename X >
constexpr X phys::units::magnitude (quantity< DX, X > const &q)
 quantity's magnitude.
 

Variables

constexpr struct phys::units::detail::magnitude_tag_t phys::units::detail::magnitude_tag
 
template<typename T >
constexpr bool phys::units::is_quantity_v = is_quantity<T>::value
 
constexpr quantity< length_d > phys::units::meter {detail::magnitude_tag, 1.0}
 
constexpr quantity< mass_d > phys::units::kilogram {detail::magnitude_tag, 1.0}
 
constexpr quantity< time_interval_d > phys::units::second {detail::magnitude_tag, 1.0}
 
constexpr quantity< electric_current_d > phys::units::ampere {detail::magnitude_tag, 1.0}
 
constexpr quantity< thermodynamic_temperature_d > phys::units::kelvin {detail::magnitude_tag, 1.0}
 
constexpr quantity< amount_of_substance_d > phys::units::mole {detail::magnitude_tag, 1.0}
 
constexpr quantity< luminous_intensity_d > phys::units::candela {detail::magnitude_tag, 1.0}
 
constexpr quantity< hepenergy_d > phys::units::electronvolt {detail::magnitude_tag, 1.0}
 
constexpr long double phys::units::quetta = 1e+30L
 
constexpr long double phys::units::ronna = 1e+27L
 
constexpr long double phys::units::yotta = 1e+24L
 
constexpr long double phys::units::zetta = 1e+21L
 
constexpr long double phys::units::exa = 1e+18L
 
constexpr long double phys::units::peta = 1e+15L
 
constexpr long double phys::units::tera = 1e+12L
 
constexpr long double phys::units::giga = 1e+9L
 
constexpr long double phys::units::mega = 1e+6L
 
constexpr long double phys::units::kilo = 1e+3L
 
constexpr long double phys::units::hecto = 1e+2L
 
constexpr long double phys::units::deka = 1e+1L
 
constexpr long double phys::units::deci = 1e-1L
 
constexpr long double phys::units::centi = 1e-2L
 
constexpr long double phys::units::milli = 1e-3L
 
constexpr long double phys::units::micro = 1e-6L
 
constexpr long double phys::units::nano = 1e-9L
 
constexpr long double phys::units::pico = 1e-12L
 
constexpr long double phys::units::femto = 1e-15L
 
constexpr long double phys::units::atto = 1e-18L
 
constexpr long double phys::units::zepto = 1e-21L
 
constexpr long double phys::units::yocto = 1e-24L
 
constexpr long double phys::units::ronto = 1e-27L
 
constexpr long double phys::units::quecto = 1e-30L
 
constexpr long double phys::units::kibi = 1024
 
constexpr long double phys::units::mebi = 1024 * kibi
 
constexpr long double phys::units::gibi = 1024 * mebi
 
constexpr long double phys::units::tebi = 1024 * gibi
 
constexpr long double phys::units::pebi = 1024 * tebi
 
constexpr long double phys::units::exbi = 1024 * pebi
 
constexpr long double phys::units::zebi = 1024 * exbi
 
constexpr long double phys::units::yobi = 1024 * zebi
 
constexpr Rep phys::units::pi {Rep(3.141592653589793238462L)}
 
constexpr Rep phys::units::percent {Rep(1) / 100}
 
constexpr quantity< mass_d > phys::units::gram {kilogram / 1000}
 
constexpr Rep phys::units::radian {Rep(1)}
 
constexpr Rep phys::units::steradian {Rep(1)}
 
constexpr quantity< force_d > phys::units::newton {meter * kilogram / square(second)}
 
constexpr quantity< pressure_d > phys::units::pascal {newton / square(meter)}
 
constexpr quantity< energy_d > phys::units::joule {newton * meter}
 
constexpr quantity< power_d > phys::units::watt {joule / second}
 
constexpr quantity< electric_charge_d > phys::units::coulomb {second * ampere}
 
constexpr quantity< electric_potential_d > phys::units::volt {watt / ampere}
 
constexpr quantity< capacitance_d > phys::units::farad {coulomb / volt}
 
constexpr quantity< electric_resistance_d > phys::units::ohm {volt / ampere}
 
constexpr quantity< electric_conductance_d > phys::units::siemens {ampere / volt}
 
constexpr quantity< magnetic_flux_d > phys::units::weber {volt * second}
 
constexpr quantity< magnetic_flux_density_d > phys::units::tesla {weber / square(meter)}
 
constexpr quantity< inductance_d > phys::units::henry {weber / ampere}
 
constexpr quantity< thermodynamic_temperature_d > phys::units::degree_celsius {kelvin}
 
constexpr quantity< luminous_flux_d > phys::units::lumen {candela * steradian}
 
constexpr quantity< illuminance_d > phys::units::lux {lumen / meter / meter}
 
constexpr quantity< activity_of_a_nuclide_d > phys::units::becquerel {1 / second}
 
constexpr quantity< absorbed_dose_d > phys::units::gray {joule / kilogram}
 
constexpr quantity< dose_equivalent_d > phys::units::sievert {joule / kilogram}
 
constexpr quantity< frequency_d > phys::units::hertz {1 / second}
 
constexpr quantity< length_d > phys::units::angstrom {Rep(1e-10L) * meter}
 
constexpr quantity< area_d > phys::units::are {Rep(1e+2L) * square(meter)}
 
constexpr quantity< pressure_d > phys::units::bar {Rep(1e+5L) * pascal}
 
constexpr quantity< area_d > phys::units::barn {Rep(1e-28L) * square(meter)}
 
constexpr quantity< activity_of_a_nuclide_d > phys::units::curie {Rep(3.7e+10L) * becquerel}
 
constexpr quantity< time_interval_d > phys::units::day {Rep(86400L) * second}
 
constexpr Rep phys::units::degree_angle {pi / 180}
 
constexpr quantity< acceleration_d > phys::units::gal {Rep(1e-2L) * meter / square(second)}
 
constexpr quantity< area_d > phys::units::hectare {Rep(1e+4L) * square(meter)}
 
constexpr quantity< time_interval_d > phys::units::hour {Rep(3600) * second}
 
constexpr quantity< speed_d > phys::units::knot {Rep(1852) / 3600 * meter / second}
 
constexpr quantity< volume_d > phys::units::liter {Rep(1e-3L) * cube(meter)}
 
constexpr quantity< time_interval_d > phys::units::minute {Rep(60) * second}
 
constexpr Rep phys::units::minute_angle {pi / 10800}
 
constexpr quantity< length_d > phys::units::mile_nautical {Rep(1852) * meter}
 
constexpr quantity< absorbed_dose_d > phys::units::rad {Rep(1e-2L) * gray}
 
constexpr quantity< dose_equivalent_d > phys::units::rem {Rep(1e-2L) * sievert}
 
constexpr quantity< exposure_d > phys::units::roentgen {Rep(2.58e-4L) * coulomb / kilogram}
 
constexpr Rep phys::units::second_angle {pi / 648000L}
 
constexpr quantity< mass_d > phys::units::ton_metric {Rep(1e+3L) * kilogram}
 
constexpr quantity< length_d > phys::units::metre {meter}
 
constexpr quantity< volume_d > phys::units::litre {liter}
 
constexpr Rep phys::units::deca {deka}
 
constexpr quantity< mass_d > phys::units::tonne {ton_metric}
 

Detailed Description

Zero-overhead dimensional analysis and unit/quantity manipulation and conversion.

Author
Michael S. Kenniston, Martin Moene
Date
7 September 2013
Since
0.4

Copyright 2013 Universiteit Leiden. All rights reserved.

Copyright (c) 2001 by Michael S. Kenniston. For the most recent version check www.xnet.com/~msk/quantity. Permission is granted to use this code without restriction so long as this copyright notice appears in all source files.

This code is provided as-is, with no warrantee of correctness.

Distributed under the Boost Software License, Version 1.0. (See accompanying file LICENSE_1_0.txt or copy at http://www.boost.org/LICENSE_1_0.txt)

Definition in file quantity.hpp.

Macro Definition Documentation

◆ QUANTITY_DEFINE_SCALING_LITERAL

#define QUANTITY_DEFINE_SCALING_LITERAL (   sfx,
  dim,
  factor 
)
Value:
constexpr quantity<dim, double> operator"" _##sfx(unsigned long long x) { \
return quantity<dim, double>(detail::magnitude_tag, factor * x); \
} \
constexpr quantity<dim, double> operator"" _##sfx(long double x) { \
return quantity<dim, double>(detail::magnitude_tag, factor * x); \
}

Definition at line 936 of file quantity.hpp.

◆ QUANTITY_DEFINE_SCALING_LITERALS

#define QUANTITY_DEFINE_SCALING_LITERALS (   pfx,
  dim,
  fact 
)
Value:
QUANTITY_DEFINE_SCALING_LITERAL(Q##pfx, dim, fact* quetta) \
QUANTITY_DEFINE_SCALING_LITERAL(R##pfx, dim, fact* ronna) \
QUANTITY_DEFINE_SCALING_LITERAL(Y##pfx, dim, fact* yotta) \
QUANTITY_DEFINE_SCALING_LITERAL(Z##pfx, dim, fact* zetta) \
QUANTITY_DEFINE_SCALING_LITERAL(E##pfx, dim, fact* exa) \
QUANTITY_DEFINE_SCALING_LITERAL(P##pfx, dim, fact* peta) \
QUANTITY_DEFINE_SCALING_LITERAL(T##pfx, dim, fact* tera) \
QUANTITY_DEFINE_SCALING_LITERAL(G##pfx, dim, fact* giga) \
QUANTITY_DEFINE_SCALING_LITERAL(M##pfx, dim, fact* mega) \
QUANTITY_DEFINE_SCALING_LITERAL(k##pfx, dim, fact* kilo) \
QUANTITY_DEFINE_SCALING_LITERAL(h##pfx, dim, fact* hecto) \
QUANTITY_DEFINE_SCALING_LITERAL(da##pfx, dim, fact* deka) \
QUANTITY_DEFINE_SCALING_LITERAL(pfx, dim, fact * 1) \
QUANTITY_DEFINE_SCALING_LITERAL(d##pfx, dim, fact* deci) \
QUANTITY_DEFINE_SCALING_LITERAL(c##pfx, dim, fact* centi) \
QUANTITY_DEFINE_SCALING_LITERAL(m##pfx, dim, fact* milli) \
QUANTITY_DEFINE_SCALING_LITERAL(u##pfx, dim, fact* micro) \
QUANTITY_DEFINE_SCALING_LITERAL(n##pfx, dim, fact* nano) \
QUANTITY_DEFINE_SCALING_LITERAL(p##pfx, dim, fact* pico) \
QUANTITY_DEFINE_SCALING_LITERAL(f##pfx, dim, fact* femto) \
QUANTITY_DEFINE_SCALING_LITERAL(a##pfx, dim, fact* atto) \
QUANTITY_DEFINE_SCALING_LITERAL(z##pfx, dim, fact* zepto) \
QUANTITY_DEFINE_SCALING_LITERAL(y##pfx, dim, fact* yocto) \
QUANTITY_DEFINE_SCALING_LITERAL(r##pfx, dim, fact* ronto) \
QUANTITY_DEFINE_SCALING_LITERAL(q##pfx, dim, fact* quecto)

Definition at line 944 of file quantity.hpp.