Skip to content

Conversion functions from degrees to radians and vice-versa. #845

New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Merged
merged 16 commits into from
Jul 9, 2024
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
64 changes: 64 additions & 0 deletions doc/specs/stdlib_math.md
Original file line number Diff line number Diff line change
Expand Up @@ -382,6 +382,70 @@ Notes: Although the angle of the complex number `0` is undefined, `argpi((0,0))`
{!example/math/example_math_argpi.f90!}
```

### `deg2rad`

#### Status

Experimental

#### Class

Elemenal function.

### Description

`deg2rad` converts phase angles from degrees to radians.

#### Syntax

`result = ` [[stdlib_math(module):deg2rad(interface)]] `(theta)`

#### Arguments

`theta`: Shall be a `real` scalar/array.

#### Return value

Returns the `real` phase angle in radians.

#### Example

```fortran
{!example/math/example_math_deg2rad.f90!}
```

### `rad2deg`

#### Status

Experimental

#### Class

Elemenal function.

### Description

`rad2deg` converts phase angles from radians to degrees.

#### Syntax

`result = ` [[stdlib_math(module):rad2deg(interface)]] `(theta)`

#### Arguments

`theta`: Shall be a `real` scalar/array.

#### Return value

Returns the `real` phase angle in degrees.

#### Example

```fortran
{!example/math/example_math_rad2deg.f90!}
```

### `is_close` function

#### Description
Expand Down
2 changes: 2 additions & 0 deletions example/math/CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -12,5 +12,7 @@ ADD_EXAMPLE(math_arange)
ADD_EXAMPLE(math_argd)
ADD_EXAMPLE(math_arg)
ADD_EXAMPLE(math_argpi)
ADD_EXAMPLE(math_deg2rad)
ADD_EXAMPLE(math_rad2deg)
ADD_EXAMPLE(math_is_close)
ADD_EXAMPLE(meshgrid)
8 changes: 8 additions & 0 deletions example/math/example_math_deg2rad.f90
Original file line number Diff line number Diff line change
@@ -0,0 +1,8 @@
program example_math_deg2rad
use stdlib_math, only: deg2rad
implicit none
print *, deg2rad(0.0) ! 0.0
print *, deg2rad(90.0) ! 1.57508
print *, deg2rad(-180.0) ! -3.1416

end program example_math_deg2rad
9 changes: 9 additions & 0 deletions example/math/example_math_rad2deg.f90
Original file line number Diff line number Diff line change
@@ -0,0 +1,9 @@
program example_math_rad2deg
use stdlib_math, only: rad2deg
use stdlib_constants, only: PI_sp
implicit none
print *, rad2deg(0.0) ! 0.0
print *, rad2deg(PI_sp / 2.0) ! 90.0
print *, rad2deg(-PI_sp) ! -3.1416

end program example_math_rad2deg
38 changes: 37 additions & 1 deletion src/stdlib_math.fypp
Original file line number Diff line number Diff line change
Expand Up @@ -15,7 +15,7 @@ module stdlib_math
#:endif
public :: DEFAULT_LINSPACE_LENGTH, DEFAULT_LOGSPACE_BASE, DEFAULT_LOGSPACE_LENGTH
public :: stdlib_meshgrid_ij, stdlib_meshgrid_xy
public :: arange, arg, argd, argpi, is_close, all_close, diff, meshgrid
public :: arange, arg, argd, argpi, deg2rad, rad2deg, is_close, all_close, diff, meshgrid

integer, parameter :: DEFAULT_LINSPACE_LENGTH = 100
integer, parameter :: DEFAULT_LOGSPACE_LENGTH = 50
Expand Down Expand Up @@ -332,6 +332,26 @@ module stdlib_math
procedure :: argpi_${k1}$
#:endfor
end interface argpi

!> Version: experimental
!>
!> `deg2rad` converts phase angles from degrees to radians.
!> ([Specification](../page/specs/stdlib_math.html#deg2rad-function))
interface deg2rad
#:for k1 in REAL_KINDS
procedure :: deg2rad_${k1}$
#:endfor
end interface deg2rad

!> Version: experimental
!>
!> `rad2deg` converts phase angles from radians to degrees.
!> ([Specification](../page/specs/stdlib_math.html#rad2deg-function))
interface rad2deg
#:for k1 in REAL_KINDS
procedure :: rad2deg_${k1}$
#:endfor
end interface rad2deg

!> Returns a boolean scalar/array where two scalar/arrays are element-wise equal within a tolerance.
!> ([Specification](../page/specs/stdlib_math.html#is_close-function))
Expand Down Expand Up @@ -453,6 +473,22 @@ contains
end function argpi_${k1}$
#:endfor

#:for k1, t1 in REAL_KINDS_TYPES
elemental function deg2rad_${k1}$(theta) result(result)
${t1}$, intent(in) :: theta
${t1}$ :: result
result = theta * PI_${k1}$ / 180.0_${k1}$

end function deg2rad_${k1}$

elemental function rad2deg_${k1}$(theta) result(result)
${t1}$, intent(in) :: theta
${t1}$ :: result
result = theta * 180.0_${k1}$ / PI_${k1}$

end function rad2deg_${k1}$
#:endfor

#:for k1, t1 in INT_KINDS_TYPES
!> Returns the greatest common divisor of two integers of kind ${k1}$
!> using the Euclidean algorithm.
Expand Down
28 changes: 27 additions & 1 deletion test/math/test_stdlib_math.fypp
Original file line number Diff line number Diff line change
Expand Up @@ -5,7 +5,7 @@
module test_stdlib_math
use testdrive, only : new_unittest, unittest_type, error_type, check, skip_test
use stdlib_math, only: clip, arg, argd, argpi, arange, is_close, all_close, diff, &
arange
arange, deg2rad, rad2deg
use stdlib_kinds, only: int8, int16, int32, int64, sp, dp, xdp, qp
implicit none

Expand Down Expand Up @@ -45,6 +45,12 @@ contains
, new_unittest("argpi-cmplx-${k1}$", test_argpi_${k1}$) &
#:endfor

!> Tests for deg2rad/rad2deg
#:for k1 in REAL_KINDS
, new_unittest("deg2rad-real-${k1}$", test_deg2rad_${k1}$) &
, new_unittest("rad2deg-real-${k1}$", test_rad2deg_${k1}$) &
#:endfor

!> Tests for `is_close` and `all_close`
#:for k1 in REAL_KINDS
, new_unittest("is_close-real-${k1}$", test_is_close_real_${k1}$) &
Expand Down Expand Up @@ -301,6 +307,26 @@ contains

end subroutine test_argpi_${k1}$
#:endfor

#:for k1 in REAL_KINDS
subroutine test_deg2rad_${k1}$(error)
type(error_type), allocatable, intent(out) :: error
real(${k1}$), parameter :: tol = sqrt(epsilon(1.0_${k1}$))

call check(error, PI_${k1}$, deg2rad(180.0_${k1}$), thr=tol)
if (allocated(error)) return

end subroutine test_deg2rad_${k1}$

subroutine test_rad2deg_${k1}$(error)
type(error_type), allocatable, intent(out) :: error
real(${k1}$), parameter :: tol = sqrt(epsilon(1.0_${k1}$))

call check(error, 180.0_${k1}$, rad2deg(PI_${k1}$), thr=tol)
if (allocated(error)) return

end subroutine test_rad2deg_${k1}$
#:endfor

#:for k1 in REAL_KINDS
subroutine test_is_close_real_${k1}$(error)
Expand Down
Loading