-
git-svn-id: https://lxsd.femto-st.fr/svn/fvn@59 b657c933-2333-4658-acf2-d3c7c2708721
test_operators.f90
2.27 KB
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
program test_matinv
use fvn_linear
implicit none
integer, parameter :: n=3
complex(kind=dp_kind),dimension(n,n) :: m1,m2,res
real(kind=dp_kind),dimension(n,n) :: rtmp,itmp
character(len=80) :: fmcmplx
integer :: i
fmcmplx='(3("(",f8.5,",",f8.5,") "))'
! initialize pseudo random generator
!call init_random_seed()
! fill real and imaginary part
call random_number(rtmp)
call random_number(itmp)
! create the complex matrix (fvn_i is defined in the fvn module)
m1=rtmp+fvn_i*itmp
write(*,*) "Matrix M1"
do i=1,n
write(*,fmcmplx) m1(i,:)
end do
call random_number(rtmp)
call random_number(itmp)
m2=rtmp+fvn_i*itmp
write(*,*)
write(*,*) "Matrix M2"
do i=1,n
write(*,fmcmplx) m2(i,:)
end do
write(*,*)
write(*,*) "M1.x.M2"
res=m1.x.m2
call write_res()
write(*,*)
write(*,*) "M1.x.M2 standard"
res=matmul(m1,m2)
call write_res()
write(*,*)
write(*,*) ".i.M1"
res=.i.m1
call write_res()
write(*,*)
write(*,*) ".i.M1 standard"
call fvn_matinv(3,m1,res)
call write_res()
write(*,*)
write(*,*) "M1.ix.M2"
res=m1.ix.m2
call write_res()
write(*,*)
write(*,*) "M1.ix.M2 standard"
call fvn_matinv(3,m1,res)
res=matmul(res,m2)
call write_res()
write(*,*)
write(*,*) "M1.xi.M2"
res=m1.xi.m2
call write_res()
write(*,*)
write(*,*) "M1.xi.M2 standard"
res=m1.xi.m2
call fvn_matinv(3,m2,res)
res=matmul(m1,res)
call write_res()
write(*,*)
write(*,*) ".t.M1"
res=.t.m1
call write_res()
write(*,*)
write(*,*) ".t.M1 standard"
res=transpose(m1)
call write_res()
write(*,*)
write(*,*) "M1.tx.M2"
res=m1.tx.m2
call write_res()
write(*,*)
write(*,*) "M1.tx.M2 standard"
res=matmul(transpose(m1),m2)
call write_res()
write(*,*)
write(*,*) "M1.xt.M2"
res=m1.xt.m2
call write_res()
write(*,*)
write(*,*) "M1.xt.M2 standard"
res=matmul(m1,transpose(m2))
call write_res()
write(*,*)
write(*,*) ".h.M1"
res=.h.m1
call write_res()
write(*,*)
write(*,*) ".h.M1 standard"
res=transpose(conjg(m1))
call write_res()
write(*,*)
write(*,*) "M1.hx.M2"
res=m1.hx.m2
call write_res()
write(*,*)
write(*,*) "M1.hx.M2 standard"
res=matmul(transpose(conjg(m1)),m2)
call write_res()
write(*,*)
write(*,*) "M1.xh.M2"
res=m1.xh.m2
call write_res()
write(*,*)
write(*,*) "M1.xh.M2 standard"
res=matmul(m1,transpose(conjg(m2)))
call write_res()
contains
subroutine write_res()
do i=1,n
write(*,fmcmplx) res(i,:)
end do
end subroutine
end program