The stdlib_constants module provides mathematical constants and the most common physical constants.
Warning: The names of the most common physical constants are kept short as they are inside a dedicated module. Nonetheless, in case of overlapping names, they can always be renamed as following:
use stdlib_constants, only: clight => c
The stdlib_codata module defines all codata (physical) constants as derived type. The module is automatically generated with a simple parser written in Python The latest codata constants were released in 2022 by the NIST All values for the codata constants are provided as double precision reals. The names are quite long and can be aliased with shorter names.
The derived type codata_constant_type defines:
4 members:
name
(string)value
(double precision real)uncertainty
(double precision real)unit
(string)2 type-bound procedures:
print
: to print the values of the constant members;to_real
: to get the value or the uncertainty to the desired precision.A module level interface to_real is available for getting the constant value or uncertainty of a constant.
to_real
- Get the constant value or its uncertainty.Experimental
Convert a codata_constant_type to a real
(at least sp
, or dp
) scalar.
Warning: Some constants cannot be converted to single precision sp
reals due to the value of the exponents.
r =
to_real (c, mold [, uncertainty])
c
: argument has intent(in)
and shall be of type codata_constant_type.
mold
: argument has intent(in)
and shall be of real
type.
Note: The type of the mold
argument defines the type of the result.
uncertainty
(optional): argument has intent(in)
and shall be of logical
type.
It specifies if the uncertainty needs to be returned instead of the value. Default to .false.
.
Returns a scalar of real
type which is either the value or the uncertainty of a codata constant.
program example_constants
use stdlib_constants, only: c, pi=>PI_dp
use stdlib_codata, only: alpha=>ALPHA_PARTICLE_ELECTRON_MASS_RATIO
use stdlib_codata_type, only : to_real
use stdlib_kinds, only: dp, sp
! Use most common physical constants defined as double precision reals
print *, "speed of light in vacuum= ", c
! Use of mathematical constants such as PI
print *, "PI as double precision real= ", pi
! Use codata_constant type for evaluating the value to the desired precision
print *, "Value of alpha... evaluated to double precision=", alpha%to_real(1.0_dp)
print *, "Uncertainty of alpha... evaluated to double precision=", alpha%to_real(1.0_sp, .true.)
print *, "Value of alpha... evaluated to single precision=", alpha%to_real(1.0_sp)
! Convert a codata constant to a real
print *, "Value of the alpha... evaluated to double precision=", to_real(alpha, 1.0_dp)
! Print out codata constant attributes: name, value, uncertainty and unit
call alpha%print()
end program example_constants