Skip to content

Commit cee99f3

Browse files
committed
Moved diag to submodule
1 parent f5538d6 commit cee99f3

4 files changed

+109
-75
lines changed

src/CMakeLists.txt

+1
Original file line numberDiff line numberDiff line change
@@ -4,6 +4,7 @@
44
set(fppFiles
55
stdlib_experimental_io.fypp
66
stdlib_experimental_linalg.fypp
7+
stdlib_experimental_linalg_diag.fypp
78
stdlib_experimental_optval.fypp
89
stdlib_experimental_stats.fypp
910
stdlib_experimental_stats_mean.fypp

src/Makefile.manual

+3-1
Original file line numberDiff line numberDiff line change
@@ -3,6 +3,7 @@ SRC = f18estop.f90 \
33
stdlib_experimental_error.f90 \
44
stdlib_experimental_io.f90 \
55
stdlib_experimental_linalg.f90 \
6+
stdlib_experimental_linalg_diag.f90 \
67
stdlib_experimental_kinds.f90 \
78
stdlib_experimental_optval.f90 \
89
stdlib_experimental_quadrature.f90 \
@@ -43,7 +44,7 @@ stdlib_experimental_io.o: \
4344
stdlib_experimental_error.o \
4445
stdlib_experimental_optval.o \
4546
stdlib_experimental_kinds.o
46-
stdlib_experimental_linalg.o: stdlib_experimental_kinds.o
47+
stdlib_experimental_linalg_diag.o: stdlib_experimental_kinds.o
4748
stdlib_experimental_optval.o: stdlib_experimental_kinds.o
4849
stdlib_experimental_quadrature.o: stdlib_experimental_kinds.o
4950
stdlib_experimental_stats_mean.o: \
@@ -62,6 +63,7 @@ stdlib_experimental_stats_var.o: \
6263
# Fortran sources that are built from fypp templates
6364
stdlib_experimental_io.f90: stdlib_experimental_io.fypp
6465
stdlib_experimental_linalg.f90: stdlib_experimental_linalg.fypp
66+
stdlib_experimental_linalg_diag.f90: stdlib_experimental_linalg_diag.fypp
6567
stdlib_experimental_quadrature.f90: stdlib_experimental_quadrature.fypp
6668
stdlib_experimental_stats.f90: stdlib_experimental_stats.fypp
6769
stdlib_experimental_stats_mean.f90: stdlib_experimental_stats_mean.fypp

src/stdlib_experimental_linalg.fypp

+25-74
Original file line numberDiff line numberDiff line change
@@ -11,19 +11,40 @@ module stdlib_experimental_linalg
1111
public :: trace
1212

1313
interface diag
14+
!
1415
! Vector to matrix
16+
!
1517
#:for k1, t1 in RCI_KINDS_TYPES
16-
module procedure diag_${t1[0]}$${k1}$
18+
module function diag_${t1[0]}$${k1}$(v) result(res)
19+
${t1}$, intent(in) :: v(:)
20+
${t1}$ :: res(size(v),size(v))
21+
end function diag_${t1[0]}$${k1}$
1722
#:endfor
23+
1824
#:for k1, t1 in RCI_KINDS_TYPES
19-
module procedure diag_${t1[0]}$${k1}$_k
25+
module function diag_${t1[0]}$${k1}$_k(v,k) result(res)
26+
${t1}$, intent(in) :: v(:)
27+
integer, intent(in) :: k
28+
${t1}$ :: res(size(v)+abs(k),size(v)+abs(k))
29+
end function diag_${t1[0]}$${k1}$_k
2030
#:endfor
31+
32+
!
2133
! Matrix to vector
34+
!
2235
#:for k1, t1 in RCI_KINDS_TYPES
23-
module procedure diag_${t1[0]}$${k1}$_mat
36+
module function diag_${t1[0]}$${k1}$_mat(A) result(res)
37+
${t1}$, intent(in) :: A(:,:)
38+
${t1}$ :: res(minval(shape(A)))
39+
end function diag_${t1[0]}$${k1}$_mat
2440
#:endfor
41+
2542
#:for k1, t1 in RCI_KINDS_TYPES
26-
module procedure diag_${t1[0]}$${k1}$_mat_k
43+
module function diag_${t1[0]}$${k1}$_mat_k(A,k) result(res)
44+
${t1}$, intent(in) :: A(:,:)
45+
integer, intent(in) :: k
46+
${t1}$ :: res(minval(shape(A))-abs(k))
47+
end function diag_${t1[0]}$${k1}$_mat_k
2748
#:endfor
2849
end interface
2950

@@ -47,76 +68,6 @@ contains
4768
end do
4869
end function eye
4970

50-
#:for k1, t1 in RCI_KINDS_TYPES
51-
function diag_${t1[0]}$${k1}$(v) result(res)
52-
${t1}$, intent(in) :: v(:)
53-
${t1}$ :: res(size(v),size(v))
54-
integer :: i
55-
res = 0
56-
do i = 1, size(v)
57-
res(i,i) = v(i)
58-
end do
59-
end function diag_${t1[0]}$${k1}$
60-
#:endfor
61-
62-
63-
#:for k1, t1 in RCI_KINDS_TYPES
64-
function diag_${t1[0]}$${k1}$_k(v,k) result(res)
65-
${t1}$, intent(in) :: v(:)
66-
integer, intent(in) :: k
67-
${t1}$ :: res(size(v)+abs(k),size(v)+abs(k))
68-
integer :: i, sz
69-
sz = size(v)
70-
res = 0
71-
if (k > 0) then
72-
do i = 1, sz
73-
res(i,k+i) = v(i)
74-
end do
75-
else if (k < 0) then
76-
do i = 1, sz
77-
res(i+abs(k),i) = v(i)
78-
end do
79-
else
80-
do i = 1, sz
81-
res(i,i) = v(i)
82-
end do
83-
end if
84-
end function diag_${t1[0]}$${k1}$_k
85-
#:endfor
86-
87-
#:for k1, t1 in RCI_KINDS_TYPES
88-
function diag_${t1[0]}$${k1}$_mat(A) result(res)
89-
${t1}$, intent(in) :: A(:,:)
90-
${t1}$ :: res(minval(shape(A)))
91-
integer :: i
92-
do i = 1, minval(shape(A))
93-
res(i) = A(i,i)
94-
end do
95-
end function diag_${t1[0]}$${k1}$_mat
96-
#:endfor
97-
98-
#:for k1, t1 in RCI_KINDS_TYPES
99-
function diag_${t1[0]}$${k1}$_mat_k(A,k) result(res)
100-
${t1}$, intent(in) :: A(:,:)
101-
integer, intent(in) :: k
102-
${t1}$ :: res(minval(shape(A))-abs(k))
103-
integer :: i, sz
104-
sz = minval(shape(A))-abs(k)
105-
if (k > 0) then
106-
do i = 1, sz
107-
res(i) = A(i,k+i)
108-
end do
109-
else if (k < 0) then
110-
do i = 1, sz
111-
res(i) = A(i+abs(k),i)
112-
end do
113-
else
114-
do i = 1, sz
115-
res(i) = A(i,i)
116-
end do
117-
end if
118-
end function diag_${t1[0]}$${k1}$_mat_k
119-
#:endfor
12071

12172
#:for k1, t1 in RCI_KINDS_TYPES
12273
function trace_${t1[0]}$${k1}$(A) result(res)
+80
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,80 @@
1+
#:include "common.fypp"
2+
#:set RCI_KINDS_TYPES = REAL_KINDS_TYPES + CMPLX_KINDS_TYPES + INT_KINDS_TYPES
3+
submodule (stdlib_experimental_linalg) stdlib_experimental_linalg_diag
4+
5+
implicit none
6+
7+
contains
8+
9+
#:for k1, t1 in RCI_KINDS_TYPES
10+
function diag_${t1[0]}$${k1}$(v) result(res)
11+
${t1}$, intent(in) :: v(:)
12+
${t1}$ :: res(size(v),size(v))
13+
integer :: i
14+
res = 0
15+
do i = 1, size(v)
16+
res(i,i) = v(i)
17+
end do
18+
end function diag_${t1[0]}$${k1}$
19+
#:endfor
20+
21+
22+
#:for k1, t1 in RCI_KINDS_TYPES
23+
function diag_${t1[0]}$${k1}$_k(v,k) result(res)
24+
${t1}$, intent(in) :: v(:)
25+
integer, intent(in) :: k
26+
${t1}$ :: res(size(v)+abs(k),size(v)+abs(k))
27+
integer :: i, sz
28+
sz = size(v)
29+
res = 0
30+
if (k > 0) then
31+
do i = 1, sz
32+
res(i,k+i) = v(i)
33+
end do
34+
else if (k < 0) then
35+
do i = 1, sz
36+
res(i+abs(k),i) = v(i)
37+
end do
38+
else
39+
do i = 1, sz
40+
res(i,i) = v(i)
41+
end do
42+
end if
43+
end function diag_${t1[0]}$${k1}$_k
44+
#:endfor
45+
46+
#:for k1, t1 in RCI_KINDS_TYPES
47+
function diag_${t1[0]}$${k1}$_mat(A) result(res)
48+
${t1}$, intent(in) :: A(:,:)
49+
${t1}$ :: res(minval(shape(A)))
50+
integer :: i
51+
do i = 1, minval(shape(A))
52+
res(i) = A(i,i)
53+
end do
54+
end function diag_${t1[0]}$${k1}$_mat
55+
#:endfor
56+
57+
#:for k1, t1 in RCI_KINDS_TYPES
58+
function diag_${t1[0]}$${k1}$_mat_k(A,k) result(res)
59+
${t1}$, intent(in) :: A(:,:)
60+
integer, intent(in) :: k
61+
${t1}$ :: res(minval(shape(A))-abs(k))
62+
integer :: i, sz
63+
sz = minval(shape(A))-abs(k)
64+
if (k > 0) then
65+
do i = 1, sz
66+
res(i) = A(i,k+i)
67+
end do
68+
else if (k < 0) then
69+
do i = 1, sz
70+
res(i) = A(i+abs(k),i)
71+
end do
72+
else
73+
do i = 1, sz
74+
res(i) = A(i,i)
75+
end do
76+
end if
77+
end function diag_${t1[0]}$${k1}$_mat_k
78+
#:endfor
79+
80+
end submodule

0 commit comments

Comments
 (0)