Skip to content

Navigation Menu

Sign in
Appearance settings

Search code, repositories, users, issues, pull requests...

Provide feedback

We read every piece of feedback, and take your input very seriously.

Saved searches

Use saved searches to filter your results more quickly

Sign up
Appearance settings

Commit 01b3fb9

Browse files
authored
[stdlib_math] add arg/argd/argpi (#498)
1 parent bae6be5 commit 01b3fb9

File tree

5 files changed

+262
-5
lines changed

5 files changed

+262
-5
lines changed

‎CHANGELOG.md

Lines changed: 5 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -16,6 +16,11 @@ Features available from the latest git source
1616
- new module `stdlib_io_npy`
1717
[#581](https://github.com/fortran-lang/stdlib/pull/581)
1818
- new procedures `save_npy`, `load_npy`
19+
- update module `stdlib_math`
20+
- new procedures `is_close` and `all_close`
21+
[#488](https://github.com/fortran-lang/stdlib/pull/488)
22+
- new procedures `arg`, `argd` and `argpi`
23+
[#498](https://github.com/fortran-lang/stdlib/pull/498)
1924

2025
Changes to existing modules
2126

‎doc/specs/stdlib_math.md

Lines changed: 123 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -389,6 +389,129 @@ program demo_math_arange
389389
end program demo_math_arange
390390
```
391391

392+
### `arg` - Computes the phase angle in radian of a complex scalar
393+
394+
#### Status
395+
396+
Experimental
397+
398+
#### Class
399+
400+
Elemental function.
401+
402+
#### Description
403+
404+
`arg` computes the phase angle (radian version) of `complex` scalar in the interval (-π,π].
405+
The angles in `θ` are such that `z = abs(z)*exp((0.0, θ))`.
406+
407+
#### Syntax
408+
409+
`result = [[stdlib_math(module):arg(interface)]](z)`
410+
411+
#### Arguments
412+
413+
`z`: Shall be a `complex` scalar/array.
414+
This is an `intent(in)` argument.
415+
416+
#### Return value
417+
418+
Returns the `real` type phase angle (radian version) of the `complex` argument `z`.
419+
420+
Notes: Although the angle of the complex number `0` is undefined, `arg((0,0))` returns the value `0`.
421+
422+
#### Example
423+
424+
```fortran
425+
program demo_math_arg
426+
use stdlib_math, only: arg
427+
print *, arg((0.0, 0.0)) !! 0.0
428+
print *, arg((3.0, 4.0)) !! 0.927
429+
print *, arg(2.0*exp((0.0, 0.5))) !! 0.5
430+
end program demo_math_arg
431+
```
432+
433+
### `argd` - Computes the phase angle in degree of a complex scalar
434+
435+
#### Status
436+
437+
Experimental
438+
439+
#### Class
440+
441+
Elemental function.
442+
443+
#### Description
444+
445+
`argd` computes the phase angle (degree version) of `complex` scalar in the interval (-180.0,180.0].
446+
The angles in `θ` are such that `z = abs(z)*exp((0.0, θ*π/180.0))`.
447+
448+
#### Syntax
449+
450+
`result = [[stdlib_math(module):argd(interface)]](z)`
451+
452+
#### Arguments
453+
454+
`z`: Shall be a `complex` scalar/array.
455+
This is an `intent(in)` argument.
456+
457+
#### Return value
458+
459+
Returns the `real` type phase angle (degree version) of the `complex` argument `z`.
460+
461+
Notes: Although the angle of the complex number `0` is undefined, `argd((0,0))` returns the value `0`.
462+
463+
#### Example
464+
465+
```fortran
466+
program demo_math_argd
467+
use stdlib_math, only: argd
468+
print *, argd((0.0, 0.0)) !! 0.0
469+
print *, argd((3.0, 4.0)) !! 53.1°
470+
print *, argd(2.0*exp((0.0, 0.5))) !! 28.64°
471+
end program demo_math_argd
472+
```
473+
474+
### `argpi` - Computes the phase angle in circular of a complex scalar
475+
476+
#### Status
477+
478+
Experimental
479+
480+
#### Class
481+
482+
Elemental function.
483+
484+
#### Description
485+
486+
`argpi` computes the phase angle (IEEE circular version) of `complex` scalar in the interval (-1.0,1.0].
487+
The angles in `θ` are such that `z = abs(z)*exp((0.0, θ*π))`.
488+
489+
#### Syntax
490+
491+
`result = [[stdlib_math(module):argpi(interface)]](z)`
492+
493+
#### Arguments
494+
495+
`z`: Shall be a `complex` scalar/array.
496+
This is an `intent(in)` argument.
497+
498+
#### Return value
499+
500+
Returns the `real` type phase angle (circular version) of the `complex` argument `z`.
501+
502+
Notes: Although the angle of the complex number `0` is undefined, `argpi((0,0))` returns the value `0`.
503+
504+
#### Example
505+
506+
```fortran
507+
program demo_math_argpi
508+
use stdlib_math, only: argpi
509+
print *, argpi((0.0, 0.0)) !! 0.0
510+
print *, argpi((3.0, 4.0)) !! 0.295
511+
print *, argpi(2.0*exp((0.0, 0.5))) !! 0.159
512+
end program demo_math_argpi
513+
```
514+
392515
### `is_close`
393516

394517
#### Description
@@ -449,7 +572,6 @@ Returns a `logical` scalar/array.
449572
program demo_math_is_close
450573
451574
use stdlib_math, only: is_close
452-
use stdlib_error, only: check
453575
real :: x(2) = [1, 2], y, NAN
454576
455577
y = -3
@@ -514,7 +636,6 @@ Returns a `logical` scalar.
514636
program demo_math_all_close
515637
516638
use stdlib_math, only: all_close
517-
use stdlib_error, only: check
518639
real :: x(2) = [1, 2], y, NAN
519640
complex :: z(4, 4)
520641

‎src/stdlib_math.fypp

Lines changed: 63 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -14,7 +14,7 @@ module stdlib_math
1414
public :: EULERS_NUMBER_QP
1515
#:endif
1616
public :: DEFAULT_LINSPACE_LENGTH, DEFAULT_LOGSPACE_BASE, DEFAULT_LOGSPACE_LENGTH
17-
public :: arange, is_close, all_close
17+
public :: arange, arg, argd, argpi, is_close, all_close
1818

1919
integer, parameter :: DEFAULT_LINSPACE_LENGTH = 100
2020
integer, parameter :: DEFAULT_LOGSPACE_LENGTH = 50
@@ -27,6 +27,11 @@ module stdlib_math
2727
real(qp), parameter :: EULERS_NUMBER_QP = exp(1.0_qp)
2828
#:endif
2929

30+
!> Useful constants `PI` for `argd/argpi`
31+
#:for k1 in REAL_KINDS
32+
real(kind=${k1}$), parameter :: PI_${k1}$ = acos(-1.0_${k1}$)
33+
#:endfor
34+
3035
interface clip
3136
#:for k1, t1 in IR_KINDS_TYPES
3237
module procedure clip_${k1}$
@@ -296,6 +301,34 @@ module stdlib_math
296301

297302
!> Version: experimental
298303
!>
304+
!> `arg` computes the phase angle in the interval (-π,π].
305+
!> ([Specification](../page/specs/stdlib_math.html#arg))
306+
interface arg
307+
#:for k1 in CMPLX_KINDS
308+
procedure :: arg_${k1}$
309+
#:endfor
310+
end interface arg
311+
312+
!> Version: experimental
313+
!>
314+
!> `argd` computes the phase angle of degree version in the interval (-180.0,180.0].
315+
!> ([Specification](../page/specs/stdlib_math.html#argd))
316+
interface argd
317+
#:for k1 in CMPLX_KINDS
318+
procedure :: argd_${k1}$
319+
#:endfor
320+
end interface argd
321+
322+
!> Version: experimental
323+
!>
324+
!> `argpi` computes the phase angle of circular version in the interval (-1.0,1.0].
325+
!> ([Specification](../page/specs/stdlib_math.html#argpi))
326+
interface argpi
327+
#:for k1 in CMPLX_KINDS
328+
procedure :: argpi_${k1}$
329+
#:endfor
330+
end interface argpi
331+
299332
!> Returns a boolean scalar/array where two scalar/arrays are element-wise equal within a tolerance.
300333
!> ([Specification](../page/specs/stdlib_math.html#is_close))
301334
interface is_close
@@ -341,6 +374,34 @@ contains
341374

342375
#:endfor
343376

377+
#:for k1, t1 in CMPLX_KINDS_TYPES
378+
elemental function arg_${k1}$(z) result(result)
379+
${t1},ドル intent(in) :: z
380+
real(${k1}$) :: result
381+
382+
result = merge(0.0_${k1},ドル atan2(z%im, z%re), z == (0.0_${k1},ドル 0.0_${k1}$))
383+
384+
end function arg_${k1}$
385+
386+
elemental function argd_${k1}$(z) result(result)
387+
${t1},ドル intent(in) :: z
388+
real(${k1}$) :: result
389+
390+
result = merge(0.0_${k1},ドル atan2(z%im, z%re), z == (0.0_${k1},ドル 0.0_${k1}$)) &
391+
*180.0_${k1}$/PI_${k1}$
392+
393+
end function argd_${k1}$
394+
395+
elemental function argpi_${k1}$(z) result(result)
396+
${t1},ドル intent(in) :: z
397+
real(${k1}$) :: result
398+
399+
result = merge(0.0_${k1},ドル atan2(z%im, z%re), z == (0.0_${k1},ドル 0.0_${k1}$)) &
400+
/PI_${k1}$
401+
402+
end function argpi_${k1}$
403+
#:endfor
404+
344405
#:for k1, t1 in INT_KINDS_TYPES
345406
!> Returns the greatest common divisor of two integers of kind ${k1}$
346407
!> using the Euclidean algorithm.
@@ -361,4 +422,5 @@ contains
361422
end function gcd_${k1}$
362423

363424
#:endfor
425+
364426
end module stdlib_math

‎src/tests/math/Makefile.manual

Lines changed: 0 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -9,5 +9,4 @@ PROGS_SRC = test_linspace.f90 test_logspace.f90 \
99
$(SRCGEN): %.f90: %.fypp ../../common.fypp
1010
fypp -I../.. $(FYPPFLAGS) $< $@
1111

12-
1312
include ../Makefile.manual.test.mk

‎src/tests/math/test_stdlib_math.fypp

Lines changed: 71 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -4,11 +4,15 @@
44

55
module test_stdlib_math
66
use testdrive, only : new_unittest, unittest_type, error_type, check, skip_test
7-
use stdlib_math, only: clip, is_close, all_close
7+
use stdlib_math, only: clip, arg, argd, argpi, arange, is_close, all_close
88
use stdlib_kinds, only: int8, int16, int32, int64, sp, dp, xdp, qp
99
implicit none
1010

1111
public :: collect_stdlib_math
12+
13+
#:for k1 in REAL_KINDS
14+
real(kind=${k1}$), parameter :: PI_${k1}$ = acos(-1.0_${k1}$)
15+
#:endfor
1216

1317
contains
1418

@@ -33,6 +37,13 @@ contains
3337
new_unittest("clip-real-quad", test_clip_rqp), &
3438
new_unittest("clip-real-quad-bounds", test_clip_rqp_bounds) &
3539

40+
!> Tests for arg/argd/argpi
41+
#:for k1 in CMPLX_KINDS
42+
, new_unittest("arg-cmplx-${k1}$", test_arg_${k1}$) &
43+
, new_unittest("argd-cmplx-${k1}$", test_argd_${k1}$) &
44+
, new_unittest("argpi-cmplx-${k1}$", test_argpi_${k1}$) &
45+
#:endfor
46+
3647
!> Tests for `is_close` and `all_close`
3748
#:for k1 in REAL_KINDS
3849
, new_unittest("is_close-real-${k1}$", test_is_close_real_${k1}$) &
@@ -211,7 +222,66 @@ contains
211222
#:endif
212223

213224
end subroutine test_clip_rqp_bounds
225+
226+
#:for k1 in CMPLX_KINDS
227+
subroutine test_arg_${k1}$(error)
228+
type(error_type), allocatable, intent(out) :: error
229+
real(${k1}$), parameter :: tol = sqrt(epsilon(1.0_${k1}$))
230+
real(${k1}$), allocatable :: theta(:)
231+
232+
#! For scalar
233+
call check(error, abs(arg(2*exp((0.0_${k1},ドル 0.5_${k1}$))) - 0.5_${k1}$) < tol, &
234+
"test_nonzero_scalar")
235+
if (allocated(error)) return
236+
call check(error, abs(arg((0.0_${k1},ドル 0.0_${k1}$)) - 0.0_${k1}$) < tol, &
237+
"test_zero_scalar")
238+
239+
#! and for array (180.0° see scalar version)
240+
theta = arange(-179.0_${k1},ドル 179.0_${k1},ドル 3.58_${k1}$)
241+
call check(error, all(abs(arg(exp(cmplx(0.0_${k1},ドル theta/180*PI_${k1},ドル ${k1}$))) - theta/180*PI_${k1}$) < tol), &
242+
"test_array")
243+
244+
end subroutine test_arg_${k1}$
245+
246+
subroutine test_argd_${k1}$(error)
247+
type(error_type), allocatable, intent(out) :: error
248+
real(${k1}$), parameter :: tol = sqrt(epsilon(1.0_${k1}$))
249+
real(${k1}$), allocatable :: theta(:)
250+
251+
#! For scalar
252+
call check(error, abs(argd((-1.0_${k1},ドル 0.0_${k1}$)) - 180.0_${k1}$) < tol, &
253+
"test_nonzero_scalar")
254+
if (allocated(error)) return
255+
call check(error, abs(argd((0.0_${k1},ドル 0.0_${k1}$)) - 0.0_${k1}$) < tol, &
256+
"test_zero_scalar")
257+
258+
#! and for array (180.0° see scalar version)
259+
theta = arange(-179.0_${k1},ドル 179.0_${k1},ドル 3.58_${k1}$)
260+
call check(error, all(abs(argd(exp(cmplx(0.0_${k1},ドル theta/180*PI_${k1},ドル ${k1}$))) - theta) < tol), &
261+
"test_array")
262+
263+
end subroutine test_argd_${k1}$
214264

265+
subroutine test_argpi_${k1}$(error)
266+
type(error_type), allocatable, intent(out) :: error
267+
real(${k1}$), parameter :: tol = sqrt(epsilon(1.0_${k1}$))
268+
real(${k1}$), allocatable :: theta(:)
269+
270+
#! For scalar
271+
call check(error, abs(argpi((-1.0_${k1},ドル 0.0_${k1}$)) - 1.0_${k1}$) < tol, &
272+
"test_nonzero_scalar")
273+
if (allocated(error)) return
274+
call check(error, abs(argpi((0.0_${k1},ドル 0.0_${k1}$)) - 0.0_${k1}$) < tol, &
275+
"test_zero_scalar")
276+
277+
#! and for array (180.0° see scalar version)
278+
theta = arange(-179.0_${k1},ドル 179.0_${k1},ドル 3.58_${k1}$)
279+
call check(error, all(abs(argpi(exp(cmplx(0.0_${k1},ドル theta/180*PI_${k1},ドル ${k1}$))) - theta/180) < tol), &
280+
"test_array")
281+
282+
end subroutine test_argpi_${k1}$
283+
#:endfor
284+
215285
#:for k1 in REAL_KINDS
216286
subroutine test_is_close_real_${k1}$(error)
217287
type(error_type), allocatable, intent(out) :: error

0 commit comments

Comments
(0)

AltStyle によって変換されたページ (->オリジナル) /