-
Notifications
You must be signed in to change notification settings - Fork 1
/
MATRIX_SPARSE.f90
1609 lines (1369 loc) · 47.5 KB
/
MATRIX_SPARSE.f90
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
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
MODULE MATRIX_SPARSE
USE SCIFOR, only: str,free_unit,zero,assert_shape,zeye,eye,arange
USE AUX_FUNCS, only: show_fmt,append
implicit none
private
!SPARSE ROW OF THE SPARSE MATRIX: note this is dynamic array
type sparse_row
integer :: size
real(8),dimension(:),allocatable :: vals
integer,dimension(:),allocatable :: cols
end type sparse_row
!SPARSE MATRIX STRUCTURE
type sparse_matrix
type(sparse_row),dimension(:),allocatable :: row
integer :: Nrow=0
integer :: Ncol=0
logical :: status=.false.
contains
procedure,pass :: init => sp_init_matrix
procedure,pass :: free => sp_free_matrix
procedure,pass :: load => sp_load_matrix
procedure,pass :: dump => sp_dump_matrix
procedure,pass :: fast_insert=> sp_fast_insert_element
procedure,pass :: insert => sp_insert_element
procedure,pass :: get => sp_get_element
procedure,pass :: show => sp_show_matrix
procedure,pass :: display => sp_display_matrix
procedure,pass :: spy => sp_spy_matrix
procedure,pass :: as_matrix => sp_as_matrix
procedure,pass :: dgr => sp_dgr_matrix
procedure,pass :: t => sp_transpose_matrix
procedure,pass :: nnz => sp_nnz_matrix
procedure,pass :: dot => sp_matmul_vector
end type sparse_matrix
interface sparse_matrix
module procedure :: sp_construct_matrix
end interface sparse_matrix
interface as_sparse
module procedure :: sp_construct_matrix
end interface as_sparse
interface as_matrix
module procedure :: sp_as_matrix
end interface as_matrix
interface sparse
module procedure :: sp_construct_matrix
end interface sparse
!EQUALITY with scalar and function (A=B, A=cmplx)
interface assignment(=)
module procedure :: sp_matrix_equal_scalar
module procedure :: sp_matrix_equal_matrix
end interface assignment(=)
!ADDITION
interface operator (+)
module procedure :: sp_plus_matrix
end interface operator (+)
!SUBTRACTION
interface operator (-)
module procedure :: sp_minus_matrix
end interface operator (-)
!SCALAR PRODUCT
interface operator(*)
module procedure :: sp_left_product_matrix_i
module procedure :: sp_left_product_matrix_d
! module procedure :: sp_left_product_matrix_c
!
module procedure :: sp_right_product_matrix_i
module procedure :: sp_right_product_matrix_d
! module procedure :: sp_right_product_matrix_c
end interface operator(*)
!SCALAR DIVISION
interface operator(/)
module procedure :: sp_right_division_matrix_i
module procedure :: sp_right_division_matrix_d
! module procedure :: sp_right_division_matrix_c
end interface operator(/)
!Matrix-Matrix PRODUCT
interface operator(.m.)
module procedure :: sp_matmul_matrix
end interface operator(.m.)
!KRONECKER PRODUCT
interface operator(.x.)
module procedure :: sp_kron_matrix
end interface operator(.x.)
interface sp_kron
module procedure :: sp_restricted_kron_matrix
end interface sp_kron
!RETURN SHAPE OF THE SPARSE MATRIX
intrinsic :: shape
interface shape
module procedure :: sp_shape_matrix
end interface shape
!EXTEND TRANSPOSE TO SPARSE MATRIX
intrinsic :: transpose
interface transpose
module procedure :: sp_transpose_matrix
end interface transpose
!EXTEND TRANSPOSE TO SPARSE MATRIX
interface hconjg
module procedure :: sp_hconjg_matrix
end interface hconjg
intrinsic :: matmul
interface matmul
module procedure :: sp_matmul_matrix
module procedure :: sp_matmul_vector
end interface matmul
interface sp_filter
module procedure :: sp_filter_matrix_1
module procedure :: sp_filter_matrix_2
end interface sp_filter
public :: sparse_matrix
public :: as_sparse
public :: as_matrix
public :: sparse
public :: assignment(=)
public :: operator(+)
public :: operator(-)
public :: operator(*)
public :: operator(/)
public :: operator(.m.)
public :: operator(.x.)
public :: sp_kron
public :: shape
public :: transpose
public :: hconjg
public :: matmul
public :: sp_eye
public :: sp_filter
contains
!+------------------------------------------------------------------+
!PURPOSE: initialize the sparse matrix list
!+------------------------------------------------------------------+
elemental subroutine sp_init_matrix(sparse,Nrow,Ncol)
class(sparse_matrix),intent(inout) :: sparse
integer,intent(in) :: Nrow,Ncol
integer :: i
!
if(sparse%status)call sparse%free()
!
sparse%Nrow=Nrow
sparse%Ncol=Ncol
!
if(allocated(sparse%row))deallocate(sparse%row)
allocate(sparse%row(Nrow))
do i=1,Nrow
sparse%row(i)%size=0
allocate(sparse%row(i)%vals(0)) !empty array
allocate(sparse%row(i)%cols(0)) !empty array
end do
sparse%status=.true.
end subroutine sp_init_matrix
function sp_construct_matrix(matrix) result(self)
real(8),dimension(:,:),intent(in) :: matrix
type(sparse_matrix) :: self
call self%load(matrix)
end function sp_construct_matrix
!+------------------------------------------------------------------+
!PURPOSE: free an entire sparse matrix
!+------------------------------------------------------------------+
elemental subroutine sp_free_matrix(sparse)
class(sparse_matrix),intent(inout) :: sparse
integer :: i
!
do i=1,sparse%Nrow
sparse%row(i)%Size = 0
if(allocated(sparse%row(i)%vals))deallocate(sparse%row(i)%vals)
if(allocated(sparse%row(i)%cols))deallocate(sparse%row(i)%cols)
enddo
if(allocated(sparse%row))deallocate(sparse%row)
!
sparse%Nrow=0
sparse%Ncol=0
sparse%status=.false.
end subroutine sp_free_matrix
!+------------------------------------------------------------------+
!PURPOSE: load a dense matrix into a sparse one
!+------------------------------------------------------------------+
subroutine sp_load_matrix(sparse,matrix)
class(sparse_matrix),intent(inout) :: sparse
real(8),dimension(:,:),intent(in) :: matrix
integer :: i,j,Ndim1,Ndim2
!
call sparse%free()
Ndim1=size(matrix,1)
Ndim2=size(matrix,2)
!
call sparse%init(Ndim1,Ndim2) !set status=T
do i=1,Ndim1
do j=1,Ndim2
if(matrix(i,j)/=zero)call sp_insert_element(sparse,matrix(i,j),i,j)
enddo
enddo
end subroutine sp_load_matrix
!+------------------------------------------------------------------+
!PURPOSE: dump a sparse matrix into a regular 2dim array
!+------------------------------------------------------------------+
subroutine sp_dump_matrix(sparse,matrix)
class(sparse_matrix),intent(in) :: sparse
real(8),dimension(:,:),intent(inout) :: matrix
integer :: i,j,Ndim1,Ndim2
Ndim1=size(matrix,1)
Ndim2=size(matrix,2)
call assert_shape(matrix,[sparse%Nrow,sparse%Ncol],"sp_dump_matrix","Matrix")
Matrix = zero
do i=1,Ndim1
do j=1,sparse%row(i)%Size
matrix(i,sparse%row(i)%cols(j)) = matrix(i,sparse%row(i)%cols(j)) + sparse%row(i)%vals(j)
enddo
enddo
end subroutine sp_dump_matrix
!+------------------------------------------------------------------+
!PURPOSE: dump a sparse matrix into a regular 2dim array
!+------------------------------------------------------------------+
function sp_as_matrix(sparse) result(matrix)
class(sparse_matrix),intent(in) :: sparse
real(8),dimension(sparse%Nrow,sparse%Ncol) :: matrix
integer :: i,j
matrix = zero
if(.not.sparse%status)return
do i=1,sparse%Nrow
do j=1,sparse%row(i)%Size
matrix(i,sparse%row(i)%cols(j)) = sparse%row(i)%vals(j)
enddo
enddo
end function sp_as_matrix
!+------------------------------------------------------------------+
!PURPOSE: insert an element value at position (i,j) in the sparse matrix
!+------------------------------------------------------------------+
subroutine sp_insert_element(sparse,value,i,j)
class(sparse_matrix),intent(inout) :: sparse
real(8),intent(in) :: value
integer,intent(in) :: i,j
integer :: column,pos
logical :: iadd
!
if(.not.sparse%status)stop "sp_insert_element: sparse.status=F"
column = j
!
iadd = .false. !check if column already exist
if(any(sparse%row(i)%cols == column))then !
pos = binary_search(sparse%row(i)%cols,column) !find the position column in %cols
iadd=.true. !set Iadd to true
endif
!
if(iadd)then !this column exists so just sum it up
sparse%row(i)%vals(pos)=sparse%row(i)%vals(pos) + value !add up value to the current one in %vals
else !this column is new. increase counter and store it
call append(sparse%row(i)%vals,value)
call append(sparse%row(i)%cols,column)
sparse%row(i)%Size = sparse%row(i)%Size + 1
endif
!
if(sparse%row(i)%Size > sparse%Ncol)stop "sp_insert_element ERROR: row%Size > sparse%Ncol"
!
end subroutine sp_insert_element
!no addition
subroutine sp_fast_insert_element(sparse,value,i,j)
class(sparse_matrix),intent(inout) :: sparse
real(8),intent(in) :: value
integer,intent(in) :: i,j
integer :: column,pos
!
column = j
!
call append(sparse%row(i)%vals,value)
call append(sparse%row(i)%cols,column)
sparse%row(i)%Size = sparse%row(i)%Size + 1
!
if(sparse%row(i)%Size > sparse%Ncol)stop "sp_insert_element ERROR: row%Size > sparse%Ncol"
!
end subroutine sp_fast_insert_element
!+------------------------------------------------------------------+
!PURPOSE: get the element value at position (i,j) in the sparse matrix
!+------------------------------------------------------------------+
function sp_get_element(sparse,i,j) result(val)
class(sparse_matrix),intent(in) :: sparse
integer,intent(in) :: i,j
real(8) :: val
integer :: pos
if(.not.sparse%status)stop "sp_get_element: sparse.status=F"
val=zero
if(.not.any(sparse%row(i)%cols==j))return
pos=binary_search(sparse%row(i)%cols,j)
val=sparse%row(i)%vals(pos)
! do pos=1,sparse%row(i)%size
! if(j==sparse%row(i)%cols(pos))value=sparse%row(i)%vals(pos)
! enddo
end function sp_get_element
!+------------------------------------------------------------------+
!PURPOSE: return the number of non-zero elements
!+------------------------------------------------------------------+
function sp_nnz_matrix(sparse) result(nnz)
class(sparse_matrix),intent(in) :: sparse
integer :: nnz
integer :: i
nnz=0
if(.not.sparse%status)return
do i=1,sparse%Nrow
nnz=nnz + sparse%row(i)%size
enddo
end function sp_nnz_matrix
function sp_filter_matrix_1(A,states) result(Ak)
class(sparse_matrix), intent(in) :: A
integer,dimension(:),intent(in) :: states
type(sparse_matrix) :: Ak
integer :: i,j,istate,jstate
real(8) :: val
!
if(.not.A%status)stop "sp_filter_matrix_1: A.status=F"
!
call Ak%free()
call Ak%init(size(states),size(states))
!
do istate = 1,size(states)
i=states(istate)
do jstate=1,size(states)
j=states(jstate)
val = A%get(i,j)
if(val==0d0)cycle
call append(Ak%row(istate)%vals,val)
call append(Ak%row(istate)%cols,jstate)
Ak%row(istate)%Size = Ak%row(istate)%Size + 1
!
enddo
enddo
end function sp_filter_matrix_1
function sp_filter_matrix_2(A,Istates,Jstates) result(Ak)
class(sparse_matrix), intent(in) :: A
integer,dimension(:),intent(in) :: Istates,Jstates
type(sparse_matrix) :: Ak
integer :: i,j,istate,jstate
real(8) :: val
!
call Ak%free()
call Ak%init(size(Istates),size(Jstates))
!
do istate = 1,size(Istates)
i=Istates(istate)
do jstate=1,size(Jstates)
j=Jstates(jstate)
val = A%get(i,j)
if(val==0d0)cycle
call append(Ak%row(istate)%vals,val)
call append(Ak%row(istate)%cols,jstate)
Ak%row(istate)%Size = Ak%row(istate)%Size + 1
!
enddo
enddo
end function sp_filter_matrix_2
!##################################################################
!##################################################################
! SHOW / SPY
!##################################################################
!##################################################################
!+------------------------------------------------------------------+
!PURPOSE: show matrix
!+------------------------------------------------------------------+
subroutine sp_show_matrix(sparse,fmt,file)
class(sparse_matrix) :: sparse
character(len=*),optional :: fmt,file
character(len=12) :: fmt_
integer :: unit_
integer :: i,j,Ns,NN
character(len=64) :: format
real(8) :: val
unit_=6
if(present(file))open(free_unit(unit_),file=str(file))
fmt_=str(show_fmt);if(present(fmt))fmt_=str(fmt) !ES10.3
format='(A1,'//str(fmt_)//',1x)'
! format='(A1,'//str(fmt_)//',A1,'//str(fmt_)//',A1,1x)'
if(.not.sparse%status)return
Ns=sparse%Nrow
do i=1,sparse%Nrow
do j=1,sparse%Ncol
val = sp_get_element(sparse,i,j)
write(unit_,"("//str(sparse%Ncol)//"F5.1)",advance='no')val
! write(unit_,"("//str(sparse%Ncol)//str(format)//")",advance='no')"(",dreal(val),",",dimag(val),")"
enddo
write(unit_,*)
enddo
if(present(file))close(unit_)
end subroutine sp_show_matrix
subroutine sp_display_matrix(sparse)
class(sparse_matrix) :: sparse
integer :: unit_
integer :: i,j,Ns
character(len=20) :: fmtR,fmtI
real(8) :: val
unit_=6
fmtI='(I10)'
fmtR='(ES10.3)'
if(.not.sparse%status)return
do i=1,sparse%Nrow
Ns = size(sparse%row(i)%cols)
if(Ns==0)cycle
do j=1,Ns
write(unit_,"(A1,2I5,A1,F7.3)",advance='no')"[",i,sparse%row(i)%cols(j),"]",sparse%row(i)%vals(j)
enddo
write(unit_,"(A1)",advance='yes')""
! write(unit_,"("//str(Ns)//"(I10))",advance='yes')sparse%row(i)%cols
! write(unit_,"("//str(Ns)//"(2F10.3))",advance='yes')sparse%row(i)%vals
write(unit_,*)
enddo
end subroutine sp_display_matrix
!+------------------------------------------------------------------+
!PURPOSE: pretty print a sparse matrix on a given unit using format fmt
!+------------------------------------------------------------------+
subroutine sp_spy_matrix(sparse,header)
class(sparse_matrix) :: sparse
character ( len = * ) :: header
integer :: N1,N2
character ( len = 255 ) :: command_filename
integer :: command_unit
character ( len = 255 ) :: data_filename
integer :: data_unit
integer :: i, j
integer :: nz_num
!
! Create data file.
!
if(.not.sparse%status)return
!
N1 = sparse%Nrow
N2 = sparse%Ncol
data_filename = trim ( header ) // '_data.dat'
open (unit=free_unit(data_unit), file = data_filename, status = 'replace' )
nz_num = 0
do i=1,N1
do j=1,sparse%row(i)%size
write(data_unit,'(2x,i6,2x,i6)') sparse%row(i)%cols(j),i
nz_num = nz_num + 1
enddo
enddo
close(data_unit)
!
! Create command file.
!
command_filename = "plot_"//str(header)//'_commands.gp'
open(unit = free_unit(command_unit), file = command_filename, status = 'replace' )
write(command_unit,'(a)') '#unset key'
write(command_unit,'(a)') 'set terminal postscript eps enhanced color font "Times-Roman,16"'
write(command_unit,'(a)') 'set output "|ps2pdf -sEPSCrop - '//str(header)//".pdf"//'"'
write(command_unit,'(a)') 'set size ratio -1'
write(command_unit,'(a)') 'set xlabel "<--- J --->"'
write(command_unit,'(a)') 'set ylabel "<--- I --->"'
write(command_unit,'(a,i6,a)')'set title "',nz_num,' nonzeros for '//str(header)//'"'
write(command_unit,'(a)') 'set timestamp'
write(command_unit,'(a)' )'plot [x=1:'//str(N1)//'] [y='//str(N2)//':1] "'//&
str(data_filename)//'" w p pt 5 ps 0.4 lc rgb "red"'
close ( unit = command_unit )
return
end subroutine sp_spy_matrix
!##################################################################
!##################################################################
! SPARSE MATRIX KRON PRODUCT
!##################################################################
!##################################################################
!+------------------------------------------------------------------+
!PURPOSE: Perform simple Kroenecker product of two sparse matrices
!AxB(1+rowB*(i-1):rowB*i,1+colB*(j-1):colB*j) = A(i,j)*B(:,:)
!+------------------------------------------------------------------+
function sp_kron_matrix(A,B) result(AxB)
type(sparse_matrix), intent(in) :: A,B
type(sparse_matrix) :: AxB
integer :: i,icol,j,k,kcol,l
integer :: indx_row,indx_col
real(8) :: value
if(.not.A%status)stop "sp_kron_matrix: A.status=F"
if(.not.B%status)stop "sp_kron_matrix: B.status=F"
call AxB%free()
call AxB%init(a%Nrow*b%Nrow,a%Ncol*b%Ncol)
do indx_row = 1,A%Nrow*B%Nrow
k = mod(indx_row,B%Nrow);if(k==0)k=B%Nrow
i = (indx_row-1)/B%Nrow+1
do icol=1,A%row(i)%size
j = A%row(i)%cols(icol)
do kcol=1,B%row(k)%size
l = B%row(k)%cols(kcol)
indx_col = l + (j-1)*B%Ncol
value = A%row(i)%vals(icol)*B%row(k)%vals(kcol)
!
call append(AxB%row(indx_row)%vals,value)
call append(AxB%row(indx_row)%cols,indx_col)
AxB%row(indx_row)%Size = AxB%row(indx_row)%Size + 1
!
enddo
enddo
enddo
end function sp_kron_matrix
function sp_restricted_kron_matrix(A,B,states) result(AxB)
type(sparse_matrix), intent(in) :: A,B
integer,dimension(:),intent(in) :: states
type(sparse_matrix) :: AxB,Ap,Bp
integer :: i,icol,j,k,kcol,l,istate,jstate
integer :: indx_row,indx_col
real(8) :: val,Aval,Bval
!
call AxB%free()
call AxB%init(size(states),size(states))
!
do istate = 1,size(states)
indx_row=states(istate)
i = (indx_row-1)/B%Nrow+1
k = mod(indx_row,B%Nrow);if(k==0)k=B%Nrow
do jstate=1,size(states)
indx_col=states(jstate)
j = (indx_col-1)/B%Ncol+1
l = mod(indx_col,B%Ncol);if(l==0)l=B%Ncol
if(&
(.not.any(A%row(i)%cols==j))&
.OR. &
(.not.any(B%row(k)%cols==l)) )cycle
Aval = A%get(i,j)
Bval = B%get(k,l)
call append(AxB%row(istate)%vals,Aval*Bval)
call append(AxB%row(istate)%cols,jstate)
AxB%row(istate)%Size = AxB%row(istate)%Size + 1
enddo
enddo
end function sp_restricted_kron_matrix
function sp_matmul_matrix(A,B) result(AxB)
type(sparse_matrix), intent(in) :: A,B
type(sparse_matrix) :: Bt,AxB
integer :: i,icol,j,jcol,k
real(8) :: value
!
if(.not.A%status)stop "sp_matmul_matrix: A.status=F"
if(.not.B%status)stop "sp_matmul_matrix: B.status=F"
!
call AxB%free
call AxB%init(a%Nrow,b%Ncol)
Bt = B%t()
!
do i=1,AxB%Nrow !==A.Nrow
!
do j=1,AxB%Ncol !==B.Ncol=Bt.Nrow
value=zero
do icol=1,A%row(i)%size
k = A%row(i)%cols(icol) !we sum over k
do jcol=1,Bt%row(j)%size
if(k/=Bt%row(j)%cols(jcol))cycle !if there are no elements in B^T.row(j) at index k cycle
value = value + A%row(i)%vals(icol)*Bt%row(j)%vals(jcol)
enddo
enddo
if(value==zero)cycle
call append(AxB%row(i)%vals,value)
call append(AxB%row(i)%cols,j)
AxB%row(i)%Size = AxB%row(i)%Size + 1
enddo
enddo
call Bt%free
end function sp_matmul_matrix
function sp_matmul_vector(H,v) result(Hv)
class(sparse_matrix), intent(in) :: H
real(8),dimension(:) :: v
real(8),dimension(:),allocatable :: Hv
integer :: Nloc
real(8) :: val
integer :: i,j,jcol
!
if(.not.H%status)stop "sp_matmul_vector: H.status=F"
!
if(allocated(Hv))deallocate(Hv)
allocate(Hv(size(v)))
Hv=zero
Nloc=size(v)
do i=1,Nloc
matmul: do jcol=1, H%row(i)%Size
val = H%row(i)%vals(jcol)
j = H%row(i)%cols(jcol)
Hv(i) = Hv(i) + val*v(j)
end do matmul
end do
end function sp_matmul_vector
!##################################################################
!##################################################################
! SPARSE MATRIX BASIC ALGEBRA
!##################################################################
!##################################################################
function sp_dgr_matrix(a) result(c)
class(sparse_matrix), intent(in) :: a
type(sparse_matrix) :: c
integer :: col
real(8) :: val
integer :: i,j
if(.not.a%status)stop "sp_dgr_matrix: A.status=F"
call c%init(a%Ncol,a%Nrow) !hconjg
do i=1,a%Nrow
do j=1,a%row(i)%size
col = a%row(i)%cols(j)
val = a%row(i)%vals(j)
! call c%insert(conjg(val),col,i)
call c%insert(val,col,i)
enddo
enddo
end function sp_dgr_matrix
function sp_transpose_matrix(a) result(c)
class(sparse_matrix), intent(in) :: a
type(sparse_matrix) :: c
integer :: col
real(8) :: val
integer :: i,j
if(.not.a%status)stop "sp_transpose_matrix: A.status=F"
call c%init(a%Ncol,a%Nrow) !tranpose
do i=1,a%Nrow
do j=1,a%row(i)%size
col = a%row(i)%cols(j)
val = a%row(i)%vals(j)
call c%insert(val,col,i)
enddo
enddo
end function sp_transpose_matrix
function sp_hconjg_matrix(a) result(c)
class(sparse_matrix), intent(in) :: a
type(sparse_matrix) :: c
integer :: col
real(8) :: val
integer :: i,j
if(.not.a%status)stop "sp_hconjg_matrix: A.status=F"
call c%init(a%Ncol,a%Nrow) !tranpose
do i=1,a%Nrow
do j=1,a%row(i)%size
col = a%row(i)%cols(j)
val = a%row(i)%vals(j)
! call c%insert(conjg(val),col,i)
call c%insert((val),col,i)
enddo
enddo
end function sp_hconjg_matrix
!+------------------------------------------------------------------+
!PURPOSE: Sparse matrix equality spA = spB. Deep copy
!+------------------------------------------------------------------+
subroutine sp_matrix_equal_matrix(a,b)
type(sparse_matrix),intent(inout) :: a
type(sparse_matrix),intent(in) :: b
integer :: col
real(8) :: val
integer :: i,j
if(.not.b%status)stop "sp_matrix_equal_matrix: B.status=F"
call a%free()
call a%init(b%Nrow,b%Ncol)
do i=1,b%Nrow
do j=1,b%row(i)%size
col = b%row(i)%cols(j)
val = b%row(i)%vals(j)
call a%insert(val,i,col)
enddo
enddo
end subroutine sp_matrix_equal_matrix
!+------------------------------------------------------------------+
!PURPOSE: Sparse matrix scalar equality spA = const.
!+------------------------------------------------------------------+
subroutine sp_matrix_equal_scalar(a,c)
type(sparse_matrix),intent(inout) :: a
real(8),intent(in) :: c
integer :: i,j
if(.not.a%status)stop "sp_matrix_equal_matrix: A.status=F"
do i=1,a%Nrow
do j=1,a%row(i)%size
a%row(i)%vals(j) = c
enddo
enddo
end subroutine sp_matrix_equal_scalar
!+------------------------------------------------------------------+
!PURPOSE: Sparse matrix addition spA + spB = spC
!+------------------------------------------------------------------+
function sp_plus_matrix(a,b) result(c)
type(sparse_matrix), intent(in) :: a,b
type(sparse_matrix) :: c
integer :: col
real(8) :: val
integer :: i,j
if(.not.a%status)stop "sp_plus_matrix error: a.status=F"
if(.not.b%status)stop "sp_plus_matrix error: b.status=F"
if(a%Nrow/=b%Nrow)stop "sp_plus_matrix error: a.Nrow != b.Nrow"
if(a%Ncol/=b%Ncol)stop "sp_plus_matrix error: a.Ncol != b.Ncol"
c=a !copy a into c
do i=1,b%Nrow
do j=1,b%row(i)%size
col = b%row(i)%cols(j)
val = b%row(i)%vals(j)
call c%insert(val,i,col)
enddo
enddo
end function sp_plus_matrix
!+------------------------------------------------------------------+
!PURPOSE: Sparse matrix difference spA - spB = spC
!+------------------------------------------------------------------+
function sp_minus_matrix(a,b) result(c)
type(sparse_matrix), intent(in) :: a,b
type(sparse_matrix) :: c
integer :: col
real(8) :: val
integer :: i,j
if(.not.a%status)stop "sp_minus_matrix error: a.status=F"
if(.not.b%status)stop "sp_minus_matrix error: b.status=F"
if(a%Nrow/=b%Nrow)stop "sp_minus_matrix error: a.Nrow != b.Nrow"
if(a%Ncol/=b%Ncol)stop "sp_minus_matrix error: a.Ncol != b.Ncol"
c=a !copy a into c
do i=1,b%Nrow
do j=1,b%row(i)%size
col = b%row(i)%cols(j)
val =-b%row(i)%vals(j)
call c%insert(val,i,col)
enddo
enddo
end function sp_minus_matrix
!+------------------------------------------------------------------+
!PURPOSE: Sparse matrix left scalar product const*spA = spC.
! type[Const]=integer(4),real(8),cmplx(8)
!+------------------------------------------------------------------+
function sp_left_product_matrix_i(C,A) result(B)
integer,intent(in) :: C
type(sparse_matrix),intent(in) :: A
type(sparse_matrix) :: B
integer :: col
real(8) :: val
integer :: i,j
if(.not.A%status)stop "sp_left_product_matrix error: A.status=F"
call b%free()
call b%init(a%Nrow,a%Ncol)
do i=1,a%Nrow
do j=1,a%row(i)%size
col = a%row(i)%cols(j)
val = a%row(i)%vals(j)*C
call b%insert(val,i,col)
enddo
enddo
end function sp_left_product_matrix_i
function sp_left_product_matrix_d(C,A) result(B)
real(8),intent(in) :: C
type(sparse_matrix),intent(in) :: A
type(sparse_matrix) :: B
integer :: col
real(8) :: val
integer :: i,j
if(.not.A%status)stop "sp_left_product_matrix error: A.status=F"
call b%free()
call b%init(a%Nrow,a%Ncol)
do i=1,a%Nrow
do j=1,a%row(i)%size
col = a%row(i)%cols(j)
val = a%row(i)%vals(j)*C
call b%insert(val,i,col)
enddo
enddo
end function sp_left_product_matrix_d
! function sp_left_product_matrix_c(C,A) result(B)
! real(8),intent(in) :: C
! type(sparse_matrix),intent(in) :: A
! type(sparse_matrix) :: B
! integer :: col
! real(8) :: val
! integer :: i,j
! if(.not.A%status)stop "sp_left_product_matrix error: A.status=F"
! call b%free()
! call b%init(a%Nrow,a%Ncol)
! do i=1,a%Nrow
! do j=1,a%row(i)%size
! col = a%row(i)%cols(j)
! val = a%row(i)%vals(j)*C
! call b%insert(val,i,col)
! enddo
! enddo
! end function sp_left_product_matrix_c
!+------------------------------------------------------------------+
!PURPOSE: Sparse matrix right scalar product spA*const = spC.
! type[Const]=integer(4),real(8),cmplx(8)
!+------------------------------------------------------------------+
function sp_right_product_matrix_i(A,C) result(B)
integer,intent(in) :: C
type(sparse_matrix),intent(in) :: A
type(sparse_matrix) :: B
integer :: col
real(8) :: val
integer :: i,j
if(.not.A%status)stop "sp_right_product_matrix error: A.status=F"
call b%free()
call b%init(a%Nrow,a%Ncol)
do i=1,a%Nrow
do j=1,a%row(i)%size
col = a%row(i)%cols(j)
val = a%row(i)%vals(j)*C
call b%insert(val,i,col)
enddo
enddo
end function sp_right_product_matrix_i
function sp_right_product_matrix_d(A,C) result(B)
real(8),intent(in) :: C
type(sparse_matrix),intent(in) :: A
type(sparse_matrix) :: B
integer :: col
real(8) :: val
integer :: i,j
if(.not.A%status)stop "sp_right_product_matrix error: A.status=F"
call b%free()
call b%init(a%Nrow,a%Ncol)
do i=1,a%Nrow
do j=1,a%row(i)%size
col = a%row(i)%cols(j)
val = a%row(i)%vals(j)*C
call b%insert(val,i,col)
enddo
enddo
end function sp_right_product_matrix_d
! function sp_right_product_matrix_c(A,C) result(B)
! real(8),intent(in) :: C
! type(sparse_matrix),intent(in) :: A
! type(sparse_matrix) :: B
! integer :: col
! real(8) :: val
! integer :: i,j
! if(.not.A%status)stop "sp_right_product_matrix error: A.status=F"
! call b%free()
! call b%init(a%Nrow,a%Ncol)
! do i=1,a%Nrow
! do j=1,a%row(i)%size
! col = a%row(i)%cols(j)
! val = a%row(i)%vals(j)*C
! call b%insert(val,i,col)
! enddo
! enddo
! end function sp_right_product_matrix_c
!+------------------------------------------------------------------+
!PURPOSE: Sparse matrix right scalar division spA/const = spC.
! type[Const]=integer(4),real(8),cmplx(8)
!+------------------------------------------------------------------+
function sp_right_division_matrix_i(A,C) result(B)
integer,intent(in) :: C
type(sparse_matrix),intent(in) :: A
type(sparse_matrix) :: B
integer :: col
real(8) :: val
integer :: i,j
if(.not.A%status)stop "sp_right_division_matrix error: A.status=F"
call b%free()
call b%init(a%Nrow,a%Ncol)
do i=1,a%Nrow
do j=1,a%row(i)%size
col = a%row(i)%cols(j)
val = a%row(i)%vals(j)/C
call b%insert(val,i,col)
enddo
enddo
end function sp_right_division_matrix_i
function sp_right_division_matrix_d(A,C) result(B)
real(8),intent(in) :: C
type(sparse_matrix),intent(in) :: A
type(sparse_matrix) :: B
integer :: col
real(8) :: val
integer :: i,j
if(.not.A%status)stop "sp_right_division_matrix error: A.status=F"
call b%free()
call b%init(a%Nrow,a%Ncol)
do i=1,a%Nrow
do j=1,a%row(i)%size
col = a%row(i)%cols(j)
val = a%row(i)%vals(j)/C
call b%insert(val,i,col)
enddo
enddo
end function sp_right_division_matrix_d
! function sp_right_division_matrix_c(A,C) result(B)
! real(8),intent(in) :: C
! type(sparse_matrix),intent(in) :: A
! type(sparse_matrix) :: B
! integer :: col
! real(8) :: val
! integer :: i,j
! if(.not.A%status)stop "sp_right_division_matrix error: A.status=F"
! call b%free()
! call b%init(a%Nrow,a%Ncol)
! do i=1,a%Nrow
! do j=1,a%row(i)%size
! col = a%row(i)%cols(j)
! val = a%row(i)%vals(j)/C
! call b%insert(val,i,col)
! enddo
! enddo
! end function sp_right_division_matrix_c
!+------------------------------------------------------------------+
!PURPOSE: Return the shape of a sparse matrix
!+------------------------------------------------------------------+
function sp_shape_matrix(self) result(shape)
class(sparse_matrix),intent(in) :: self
integer,dimension(2) :: shape